Line data Source code
1 : // The libMesh Finite Element Library.
2 : // Copyright (C) 2002-2025 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 :
4 : // This library is free software; you can redistribute it and/or
5 : // modify it under the terms of the GNU Lesser General Public
6 : // License as published by the Free Software Foundation; either
7 : // version 2.1 of the License, or (at your option) any later version.
8 :
9 : // This library is distributed in the hope that it will be useful,
10 : // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 : // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 : // Lesser General Public License for more details.
13 :
14 : // You should have received a copy of the GNU Lesser General Public
15 : // License along with this library; if not, write to the Free Software
16 : // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 :
18 :
19 :
20 : // Local includes
21 : #include "libmesh/dof_map.h"
22 : #include "libmesh/elem.h"
23 : #include "libmesh/enum_to_string.h"
24 : #include "libmesh/fe.h"
25 : #include "libmesh/fe_interface.h"
26 : #include "libmesh/fe_macro.h"
27 : #include "libmesh/tensor_value.h"
28 :
29 :
30 : namespace libMesh
31 : {
32 :
33 :
34 0 : LIBMESH_DEFAULT_VECTORIZED_FE(0,NEDELEC_ONE)
35 0 : LIBMESH_DEFAULT_VECTORIZED_FE(1,NEDELEC_ONE)
36 2185740 : LIBMESH_DEFAULT_VECTORIZED_FE(2,NEDELEC_ONE)
37 1718794 : LIBMESH_DEFAULT_VECTORIZED_FE(3,NEDELEC_ONE)
38 :
39 :
40 : // Anonymous namespace for local helper functions
41 : namespace {
42 82416 : void nedelec_one_nodal_soln(const Elem * elem,
43 : const Order order,
44 : const std::vector<Number> & elem_soln,
45 : const int dim,
46 : const int vdim,
47 : std::vector<Number> & nodal_soln,
48 : const bool add_p_level)
49 : {
50 82416 : const unsigned int n_nodes = elem->n_nodes();
51 82416 : const ElemType elem_type = elem->type();
52 :
53 89284 : const Order totalorder = order + add_p_level*elem->p_level();
54 :
55 82416 : nodal_soln.resize(n_nodes*vdim);
56 :
57 6868 : FEType p_refined_fe_type(totalorder, NEDELEC_ONE);
58 :
59 82416 : if (elem_type != TRI6 && elem_type != TRI7 &&
60 43168 : elem_type != QUAD8 && elem_type != QUAD9 &&
61 18432 : elem_type != TET10 && elem_type != TET14 &&
62 1472 : elem_type != HEX20 && elem_type != HEX27)
63 0 : libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(elem_type) << " selected for NEDELEC_ONE FE family!");
64 :
65 82416 : const unsigned int n_sf = FEInterface::n_shape_functions(p_refined_fe_type, elem, false);
66 :
67 13736 : std::vector<Point> refspace_nodes;
68 82416 : FEVectorBase::get_refspace_nodes(elem_type,refspace_nodes);
69 6868 : libmesh_assert_equal_to (refspace_nodes.size(), n_nodes);
70 6868 : libmesh_assert_equal_to (elem_soln.size(), n_sf);
71 :
72 : // Need to create new fe object so the shape function has the FETransformation
73 : // applied to it.
74 89284 : std::unique_ptr<FEVectorBase> vis_fe = FEVectorBase::build(dim, p_refined_fe_type);
75 :
76 6868 : const std::vector<std::vector<RealGradient>> & vis_phi = vis_fe->get_phi();
77 :
78 82416 : vis_fe->reinit(elem,&refspace_nodes);
79 :
80 : // Zero before summation
81 6868 : std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
82 :
83 759120 : for (unsigned int n = 0; n < n_nodes; n++)
84 : // u = Sum (u_i phi_i)
85 5067840 : for (unsigned int i=0; i<n_sf; i++)
86 14712480 : for (int d = 0; d < vdim; d++)
87 13761792 : nodal_soln[vdim*n+d] += elem_soln[i]*(vis_phi[i][n](d));
88 :
89 89284 : return;
90 68680 : } // nedelec_one_nodal_soln
91 :
92 :
93 :
94 10823924 : unsigned int nedelec_one_n_dofs(const ElemType t, const Order o)
95 : {
96 10330832 : libmesh_assert_greater (o, 0);
97 10823924 : switch (t)
98 : {
99 4189024 : case TRI6:
100 : case TRI7:
101 4406408 : return o*(o+2);
102 3221336 : case QUAD8:
103 : case QUAD9:
104 3337116 : return 2*o*(o+1);
105 2772608 : case TET10:
106 2772608 : libmesh_assert_less (o, 2);
107 : libmesh_fallthrough();
108 : case TET14:
109 2934656 : return o*(o+2)*(o+3)/2;
110 71244 : case HEX20:
111 71244 : libmesh_assert_less (o, 2);
112 : libmesh_fallthrough();
113 : case HEX27:
114 159568 : return 3*o*(o+1)*(o+1);
115 0 : case INVALID_ELEM:
116 0 : return 0;
117 0 : default:
118 0 : libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for NEDELEC_ONE FE family!");
119 : }
120 : }
121 :
122 :
123 :
124 569542 : unsigned int nedelec_one_n_dofs(const Elem * e, const Order o)
125 : {
126 121604 : libmesh_assert(e);
127 614696 : return nedelec_one_n_dofs(e->type(), o);
128 : }
129 :
130 :
131 :
132 9604389 : unsigned int nedelec_one_n_dofs_at_node(const ElemType t,
133 : const Order o,
134 : const unsigned int n)
135 : {
136 813070 : libmesh_assert_greater (o, 0);
137 9604389 : switch (t)
138 : {
139 3071292 : case TRI6:
140 : case TRI7:
141 : {
142 2801188 : switch (n)
143 : {
144 135036 : case 0:
145 : case 1:
146 : case 2:
147 135036 : return 0;
148 1333740 : case 3:
149 : case 4:
150 : case 5:
151 1333740 : return o;
152 19280 : case 6:
153 19280 : libmesh_assert_equal_to(t, TRI7);
154 19280 : return 0;
155 0 : default:
156 0 : libmesh_error_msg("ERROR: Invalid node ID " << n);
157 : }
158 : }
159 2725959 : case QUAD8:
160 : case QUAD9:
161 : {
162 2494001 : switch (n)
163 : {
164 116728 : case 0:
165 : case 1:
166 : case 2:
167 : case 3:
168 116728 : return 0;
169 1181444 : case 4:
170 : case 5:
171 : case 6:
172 : case 7:
173 1181444 : return o;
174 15302 : case 8:
175 15302 : libmesh_assert_equal_to(t, QUAD9);
176 15302 : return 0;
177 0 : default:
178 0 : libmesh_error_msg("ERROR: Invalid node ID " << n);
179 : }
180 : }
181 3127630 : case TET10:
182 259328 : libmesh_assert_less (o, 2);
183 : libmesh_fallthrough();
184 : case TET14:
185 : {
186 2868302 : switch (n)
187 : {
188 111104 : case 0:
189 : case 1:
190 : case 2:
191 : case 3:
192 111104 : return 0;
193 1773774 : case 4:
194 : case 5:
195 : case 6:
196 : case 7:
197 : case 8:
198 : case 9:
199 1773774 : return o;
200 0 : case 10:
201 : case 11:
202 : case 12:
203 : case 13:
204 0 : libmesh_assert_equal_to(t, TET14);
205 0 : return o*(o-1);
206 0 : default:
207 0 : libmesh_error_msg("ERROR: Invalid node ID " << n);
208 : }
209 : }
210 649972 : case HEX20:
211 22144 : libmesh_assert_less (o, 2);
212 : libmesh_fallthrough();
213 : case HEX27:
214 : {
215 627828 : switch (n)
216 : {
217 18944 : case 0:
218 : case 1:
219 : case 2:
220 : case 3:
221 : case 4:
222 : case 5:
223 : case 6:
224 : case 7:
225 18944 : return 0;
226 327288 : case 8:
227 : case 9:
228 : case 10:
229 : case 11:
230 : case 12:
231 : case 13:
232 : case 14:
233 : case 15:
234 : case 16:
235 : case 17:
236 : case 18:
237 : case 19:
238 327288 : return o;
239 6336 : case 20:
240 : case 21:
241 : case 22:
242 : case 23:
243 : case 24:
244 : case 25:
245 6336 : libmesh_assert_equal_to(t, HEX27);
246 81816 : return 2*o*(o-1);
247 1056 : case 26:
248 1056 : libmesh_assert_equal_to(t, HEX27);
249 1056 : return 0;
250 0 : default:
251 0 : libmesh_error_msg("ERROR: Invalid node ID " << n);
252 : }
253 : }
254 0 : case INVALID_ELEM:
255 0 : return 0;
256 0 : default:
257 0 : libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for NEDELEC_ONE FE family!");
258 : }
259 : }
260 :
261 :
262 :
263 1733083 : unsigned int nedelec_one_n_dofs_at_node(const Elem & e,
264 : const Order o,
265 : const unsigned int n)
266 : {
267 1898099 : return nedelec_one_n_dofs_at_node(e.type(), o, n);
268 : }
269 :
270 :
271 :
272 1062844 : unsigned int nedelec_one_n_dofs_per_elem(const ElemType t,
273 : const Order o)
274 : {
275 90394 : libmesh_assert_greater (o, 0);
276 1062844 : switch (t)
277 : {
278 38596 : case TRI6:
279 : case TRI7:
280 444580 : return o*(o-1);
281 24982 : case QUAD8:
282 : case QUAD9:
283 295361 : return 2*o*(o-1);
284 24704 : case TET10:
285 24704 : libmesh_assert_less (o, 2);
286 : libmesh_fallthrough();
287 : case TET14:
288 320333 : return o*(o-1)*(o-2)/2;
289 1056 : case HEX20:
290 1056 : libmesh_assert_less (o, 2);
291 : libmesh_fallthrough();
292 : case HEX27:
293 27274 : return 3*o*(o-1)*(o-1);
294 0 : case INVALID_ELEM:
295 0 : return 0;
296 0 : default:
297 0 : libmesh_error_msg("ERROR: Invalid ElemType " << Utility::enum_to_string(t) << " selected for NEDELEC_ONE FE family!");
298 : }
299 : }
300 :
301 :
302 :
303 972450 : unsigned int nedelec_one_n_dofs_per_elem(const Elem & e,
304 : const Order o)
305 : {
306 1062844 : return nedelec_one_n_dofs_per_elem(e.type(), o);
307 : }
308 :
309 :
310 :
311 : #ifdef LIBMESH_ENABLE_AMR
312 0 : void nedelec_one_compute_constraints (DofConstraints & /*constraints*/,
313 : DofMap & /*dof_map*/,
314 : const unsigned int /*variable_number*/,
315 : const Elem * libmesh_dbg_var(elem),
316 : const unsigned Dim)
317 : {
318 : // Only constrain elements in 2,3D.
319 0 : if (Dim == 1)
320 0 : return;
321 :
322 0 : libmesh_assert(elem);
323 :
324 0 : libmesh_not_implemented();
325 : } // nedelec_one_compute_constraints()
326 : #endif // #ifdef LIBMESH_ENABLE_AMR
327 :
328 : } // anonymous namespace
329 :
330 : #define NEDELEC_LOW_D_ERROR_MESSAGE \
331 : libmesh_error_msg("ERROR: This method makes no sense for low-D elements!");
332 :
333 :
334 : // Do full-specialization for every dimension, instead
335 : // of explicit instantiation at the end of this file.
336 : template <>
337 0 : void FE<0,NEDELEC_ONE>::nodal_soln(const Elem *,
338 : const Order,
339 : const std::vector<Number> &,
340 : std::vector<Number> &,
341 : bool,
342 : const unsigned)
343 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
344 :
345 : template <>
346 0 : void FE<1,NEDELEC_ONE>::nodal_soln(const Elem *,
347 : const Order,
348 : const std::vector<Number> &,
349 : std::vector<Number> &,
350 : bool,
351 : const unsigned)
352 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
353 :
354 : template <>
355 62448 : void FE<2,NEDELEC_ONE>::nodal_soln(const Elem * elem,
356 : const Order order,
357 : const std::vector<Number> & elem_soln,
358 : std::vector<Number> & nodal_soln,
359 : const bool add_p_level,
360 : const unsigned vdim)
361 62448 : { nedelec_one_nodal_soln(elem, order, elem_soln, 2 /*dim*/, vdim, nodal_soln, add_p_level); }
362 :
363 : template <>
364 19968 : void FE<3,NEDELEC_ONE>::nodal_soln(const Elem * elem,
365 : const Order order,
366 : const std::vector<Number> & elem_soln,
367 : std::vector<Number> & nodal_soln,
368 : const bool add_p_level,
369 : const unsigned)
370 19968 : { nedelec_one_nodal_soln(elem, order, elem_soln, 3 /*dim*/, 3 /*vdim*/, nodal_soln, add_p_level); }
371 :
372 0 : LIBMESH_FE_SIDE_NODAL_SOLN(NEDELEC_ONE)
373 :
374 :
375 : // Do full-specialization for every dimension, instead
376 : // of explicit instantiation at the end of this function.
377 0 : template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs(const ElemType, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE }
378 0 : template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs(const ElemType, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE }
379 7327116 : template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs(const ElemType t, const Order o) { return nedelec_one_n_dofs(t, o); }
380 2882112 : template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs(const ElemType t, const Order o) { return nedelec_one_n_dofs(t, o); }
381 :
382 0 : template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs(const Elem *, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE }
383 0 : template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs(const Elem *, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE }
384 416408 : template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs(const Elem * e, const Order o) { return nedelec_one_n_dofs(e, o); }
385 198288 : template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs(const Elem * e, const Order o) { return nedelec_one_n_dofs(e, o); }
386 :
387 0 : template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { NEDELEC_LOW_D_ERROR_MESSAGE }
388 0 : template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { NEDELEC_LOW_D_ERROR_MESSAGE }
389 4650728 : template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return nedelec_one_n_dofs_at_node(t, o, n); }
390 3055562 : template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return nedelec_one_n_dofs_at_node(t, o, n); }
391 :
392 0 : template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { NEDELEC_LOW_D_ERROR_MESSAGE }
393 0 : template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { NEDELEC_LOW_D_ERROR_MESSAGE }
394 1146523 : template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return nedelec_one_n_dofs_at_node(e, o, n); }
395 751576 : template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return nedelec_one_n_dofs_at_node(e, o, n); }
396 :
397 0 : template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs_per_elem(const ElemType, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE }
398 0 : template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs_per_elem(const ElemType, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE }
399 0 : template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs_per_elem(const ElemType t, const Order o) { return nedelec_one_n_dofs_per_elem(t, o); }
400 0 : template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs_per_elem(const ElemType t, const Order o) { return nedelec_one_n_dofs_per_elem(t, o); }
401 :
402 0 : template <> unsigned int FE<0,NEDELEC_ONE>::n_dofs_per_elem(const Elem &, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE }
403 0 : template <> unsigned int FE<1,NEDELEC_ONE>::n_dofs_per_elem(const Elem &, const Order) { NEDELEC_LOW_D_ERROR_MESSAGE }
404 739941 : template <> unsigned int FE<2,NEDELEC_ONE>::n_dofs_per_elem(const Elem & e, const Order o) { return nedelec_one_n_dofs_per_elem(e, o); }
405 322903 : template <> unsigned int FE<3,NEDELEC_ONE>::n_dofs_per_elem(const Elem & e, const Order o) { return nedelec_one_n_dofs_per_elem(e, o); }
406 :
407 : // Nedelec first type FEMs are always tangentially continuous
408 0 : template <> FEContinuity FE<0,NEDELEC_ONE>::get_continuity() const { NEDELEC_LOW_D_ERROR_MESSAGE }
409 0 : template <> FEContinuity FE<1,NEDELEC_ONE>::get_continuity() const { NEDELEC_LOW_D_ERROR_MESSAGE }
410 0 : template <> FEContinuity FE<2,NEDELEC_ONE>::get_continuity() const { return H_CURL; }
411 0 : template <> FEContinuity FE<3,NEDELEC_ONE>::get_continuity() const { return H_CURL; }
412 :
413 : // Nedelec first type FEMs are not hierarchic
414 0 : template <> bool FE<0,NEDELEC_ONE>::is_hierarchic() const { NEDELEC_LOW_D_ERROR_MESSAGE }
415 0 : template <> bool FE<1,NEDELEC_ONE>::is_hierarchic() const { NEDELEC_LOW_D_ERROR_MESSAGE }
416 0 : template <> bool FE<2,NEDELEC_ONE>::is_hierarchic() const { return false; }
417 0 : template <> bool FE<3,NEDELEC_ONE>::is_hierarchic() const { return false; }
418 :
419 : // Nedelec first type FEM shapes always need to be reinit'ed (because of orientation dependence)
420 0 : template <> bool FE<0,NEDELEC_ONE>::shapes_need_reinit() const { NEDELEC_LOW_D_ERROR_MESSAGE }
421 0 : template <> bool FE<1,NEDELEC_ONE>::shapes_need_reinit() const { NEDELEC_LOW_D_ERROR_MESSAGE }
422 512764 : template <> bool FE<2,NEDELEC_ONE>::shapes_need_reinit() const { return true; }
423 257338 : template <> bool FE<3,NEDELEC_ONE>::shapes_need_reinit() const { return true; }
424 :
425 : #ifdef LIBMESH_ENABLE_AMR
426 : template <>
427 0 : void FE<0,NEDELEC_ONE>::compute_constraints (DofConstraints &,
428 : DofMap &,
429 : const unsigned int,
430 : const Elem *)
431 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
432 :
433 : template <>
434 0 : void FE<1,NEDELEC_ONE>::compute_constraints (DofConstraints &,
435 : DofMap &,
436 : const unsigned int,
437 : const Elem *)
438 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
439 :
440 : template <>
441 0 : void FE<2,NEDELEC_ONE>::compute_constraints (DofConstraints & constraints,
442 : DofMap & dof_map,
443 : const unsigned int variable_number,
444 : const Elem * elem)
445 0 : { nedelec_one_compute_constraints(constraints, dof_map, variable_number, elem, /*Dim=*/2); }
446 :
447 : template <>
448 0 : void FE<3,NEDELEC_ONE>::compute_constraints (DofConstraints & constraints,
449 : DofMap & dof_map,
450 : const unsigned int variable_number,
451 : const Elem * elem)
452 0 : { nedelec_one_compute_constraints(constraints, dof_map, variable_number, elem, /*Dim=*/3); }
453 : #endif // LIBMESH_ENABLE_AMR
454 :
455 : // Specialize useless shape function methods
456 : template <>
457 0 : RealGradient FE<0,NEDELEC_ONE>::shape(const ElemType, const Order,const unsigned int,const Point &)
458 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
459 : template <>
460 0 : RealGradient FE<0,NEDELEC_ONE>::shape(const Elem *,const Order,const unsigned int,const Point &,const bool)
461 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
462 : template <>
463 0 : RealGradient FE<0,NEDELEC_ONE>::shape_deriv(const ElemType, const Order,const unsigned int,
464 : const unsigned int,const Point &)
465 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
466 : template <>
467 0 : RealGradient FE<0,NEDELEC_ONE>::shape_deriv(const Elem *,const Order,const unsigned int,
468 : const unsigned int,const Point &,const bool)
469 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
470 :
471 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
472 : template <>
473 0 : RealGradient FE<0,NEDELEC_ONE>::shape_second_deriv(const ElemType, const Order,const unsigned int,
474 : const unsigned int,const Point &)
475 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
476 : template <>
477 0 : RealGradient FE<0,NEDELEC_ONE>::shape_second_deriv(const Elem *,const Order,const unsigned int,
478 : const unsigned int,const Point &,const bool)
479 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
480 :
481 : #endif
482 :
483 : template <>
484 0 : RealGradient FE<1,NEDELEC_ONE>::shape(const ElemType, const Order,const unsigned int,const Point &)
485 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
486 : template <>
487 0 : RealGradient FE<1,NEDELEC_ONE>::shape(const Elem *,const Order,const unsigned int,const Point &,const bool)
488 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
489 : template <>
490 0 : RealGradient FE<1,NEDELEC_ONE>::shape_deriv(const ElemType, const Order,const unsigned int,
491 : const unsigned int,const Point &)
492 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
493 : template <>
494 0 : RealGradient FE<1,NEDELEC_ONE>::shape_deriv(const Elem *,const Order,const unsigned int,
495 : const unsigned int,const Point &,const bool)
496 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
497 :
498 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
499 : template <>
500 0 : RealGradient FE<1,NEDELEC_ONE>::shape_second_deriv(const ElemType, const Order,const unsigned int,
501 : const unsigned int,const Point &)
502 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
503 : template <>
504 0 : RealGradient FE<1,NEDELEC_ONE>::shape_second_deriv(const Elem *,const Order,const unsigned int,
505 : const unsigned int,const Point &,const bool)
506 0 : { NEDELEC_LOW_D_ERROR_MESSAGE }
507 :
508 : #endif
509 :
510 : } // namespace libMesh
|