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 : // Local includes
19 : #include "libmesh/dof_map.h"
20 : #include "libmesh/elem.h"
21 : #include "libmesh/fe.h"
22 : #include "libmesh/fe_interface.h"
23 : #include "libmesh/fe_macro.h"
24 : #include "libmesh/tensor_value.h"
25 :
26 :
27 : namespace libMesh
28 : {
29 :
30 :
31 0 : LIBMESH_DEFAULT_VECTORIZED_FE(0,MONOMIAL_VEC)
32 0 : LIBMESH_DEFAULT_VECTORIZED_FE(1,MONOMIAL_VEC)
33 1748970 : LIBMESH_DEFAULT_VECTORIZED_FE(2,MONOMIAL_VEC)
34 0 : LIBMESH_DEFAULT_VECTORIZED_FE(3,MONOMIAL_VEC)
35 :
36 :
37 : // ------------------------------------------------------------
38 : // Vector monomial specific implementations
39 :
40 : // Anonymous namespace for local helper functions
41 : namespace
42 : {
43 : void
44 2475 : monomial_vec_nodal_soln(const Elem * elem,
45 : const Order order,
46 : const std::vector<Number> & elem_soln,
47 : const int dim,
48 : std::vector<Number> & nodal_soln,
49 : const bool add_p_level)
50 : {
51 2475 : const unsigned int n_nodes = elem->n_nodes();
52 :
53 2475 : const ElemType elem_type = elem->type();
54 :
55 2475 : nodal_soln.resize(dim * n_nodes);
56 :
57 2700 : const Order totalorder = order + add_p_level*elem->p_level();
58 :
59 2475 : switch (totalorder)
60 : {
61 : // Constant shape functions
62 0 : case CONSTANT:
63 : {
64 0 : libmesh_assert_equal_to(elem_soln.size(), static_cast<unsigned int>(dim));
65 0 : switch (dim)
66 : {
67 0 : case 2:
68 : case 3:
69 : {
70 0 : for (unsigned int n = 0; n < n_nodes; n++)
71 0 : std::copy(elem_soln.begin(), elem_soln.end(), nodal_soln.begin() + dim*n);
72 0 : return;
73 : }
74 0 : default:
75 0 : libmesh_error_msg(
76 : "The monomial_vec_nodal_soln helper should only be called for 2 and 3 dimensions");
77 : }
78 : }
79 :
80 : // For other orders, do interpolation at the nodes
81 : // explicitly.
82 225 : default:
83 : {
84 : // FEType object to be passed to various FEInterface functions below.
85 225 : FEType fe_type(order, MONOMIAL);
86 :
87 2475 : const unsigned int n_sf = FEInterface::n_shape_functions(fe_type, elem);
88 :
89 450 : std::vector<Point> refspace_nodes;
90 2475 : FEBase::get_refspace_nodes(elem_type, refspace_nodes);
91 225 : libmesh_assert_equal_to(refspace_nodes.size(), n_nodes);
92 225 : libmesh_assert_equal_to(elem_soln.size(), n_sf * dim);
93 :
94 : // Zero before summation
95 225 : std::fill(nodal_soln.begin(), nodal_soln.end(), 0);
96 :
97 7425 : for (unsigned int d = 0; d < static_cast<unsigned int>(dim); d++)
98 24750 : for (unsigned int n = 0; n < n_nodes; n++)
99 : // u_i = Sum (alpha_i phi_i)
100 79200 : for (unsigned int i = 0; i < n_sf; i++)
101 64800 : nodal_soln[d + dim * n] += elem_soln[d + dim * i] *
102 64800 : FEInterface::shape(fe_type, elem, i, refspace_nodes[n]);
103 :
104 225 : return;
105 : } // default
106 : } // switch
107 : }
108 : } // anonymous namespace
109 :
110 : // Do full-specialization for every dimension, instead
111 : // of explicit instantiation at the end of this file.
112 : // This could be macro-ified so that it fits on one line...
113 0 : LIBMESH_FE_NODAL_SOLN_DIM(MONOMIAL_VEC, (FE<0, MONOMIAL>::nodal_soln), 0)
114 0 : LIBMESH_FE_NODAL_SOLN_DIM(MONOMIAL_VEC, (FE<1, MONOMIAL>::nodal_soln), 1)
115 :
116 : template <>
117 : void
118 2475 : FE<2, MONOMIAL_VEC>::nodal_soln(const Elem * elem,
119 : const Order order,
120 : const std::vector<Number> & elem_soln,
121 : std::vector<Number> & nodal_soln,
122 : const bool add_p_level,
123 : const unsigned)
124 : {
125 2475 : monomial_vec_nodal_soln(elem, order, elem_soln, 2 /*dim*/, nodal_soln, add_p_level);
126 2475 : }
127 :
128 : template <>
129 : void
130 0 : FE<3, MONOMIAL_VEC>::nodal_soln(const Elem * elem,
131 : const Order order,
132 : const std::vector<Number> & elem_soln,
133 : std::vector<Number> & nodal_soln,
134 : const bool add_p_level,
135 : const unsigned)
136 : {
137 0 : monomial_vec_nodal_soln(elem, order, elem_soln, 3 /*dim*/, nodal_soln, add_p_level);
138 0 : }
139 :
140 0 : LIBMESH_FE_SIDE_NODAL_SOLN(MONOMIAL_VEC)
141 :
142 :
143 : // Specialize for shape function routines by leveraging scalar MONOMIAL elements
144 :
145 : // 0-D
146 : template <>
147 : RealVectorValue
148 0 : FE<0, MONOMIAL_VEC>::shape(const ElemType type,
149 : const Order order,
150 : const unsigned int i,
151 : const Point & p)
152 : {
153 0 : Real value = FE<0, MONOMIAL>::shape(type, order, i, p);
154 0 : return libMesh::RealVectorValue(value);
155 : }
156 : template <>
157 : RealVectorValue
158 0 : FE<0, MONOMIAL_VEC>::shape_deriv(const ElemType type,
159 : const Order order,
160 : const unsigned int i,
161 : const unsigned int j,
162 : const Point & p)
163 : {
164 0 : Real value = FE<0, MONOMIAL>::shape_deriv(type, order, i, j, p);
165 0 : return libMesh::RealVectorValue(value);
166 : }
167 :
168 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
169 :
170 : template <>
171 : RealVectorValue
172 0 : FE<0, MONOMIAL_VEC>::shape_second_deriv(const ElemType type,
173 : const Order order,
174 : const unsigned int i,
175 : const unsigned int j,
176 : const Point & p)
177 : {
178 0 : Real value = FE<0, MONOMIAL>::shape_second_deriv(type, order, i, j, p);
179 0 : return libMesh::RealVectorValue(value);
180 : }
181 : #endif
182 :
183 : // 1-D
184 : template <>
185 : RealVectorValue
186 0 : FE<1, MONOMIAL_VEC>::shape(const ElemType type,
187 : const Order order,
188 : const unsigned int i,
189 : const Point & p)
190 : {
191 0 : Real value = FE<1, MONOMIAL>::shape(type, order, i, p);
192 0 : return libMesh::RealVectorValue(value);
193 : }
194 : template <>
195 : RealVectorValue
196 0 : FE<1, MONOMIAL_VEC>::shape_deriv(const ElemType type,
197 : const Order order,
198 : const unsigned int i,
199 : const unsigned int j,
200 : const Point & p)
201 : {
202 0 : Real value = FE<1, MONOMIAL>::shape_deriv(type, order, i, j, p);
203 0 : return libMesh::RealVectorValue(value);
204 : }
205 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
206 : template <>
207 : RealVectorValue
208 0 : FE<1, MONOMIAL_VEC>::shape_second_deriv(const ElemType type,
209 : const Order order,
210 : const unsigned int i,
211 : const unsigned int j,
212 : const Point & p)
213 : {
214 0 : Real value = FE<1, MONOMIAL>::shape_second_deriv(type, order, i, j, p);
215 0 : return libMesh::RealVectorValue(value);
216 : }
217 :
218 : #endif
219 :
220 : // 2-D
221 : template <>
222 : RealVectorValue
223 0 : FE<2, MONOMIAL_VEC>::shape(const ElemType type,
224 : const Order order,
225 : const unsigned int i,
226 : const Point & p)
227 : {
228 0 : Real value = FE<2, MONOMIAL>::shape(type, order, i / 2, p);
229 :
230 0 : switch (i % 2)
231 : {
232 0 : case 0:
233 0 : return libMesh::RealVectorValue(value);
234 :
235 0 : case 1:
236 0 : return libMesh::RealVectorValue(Real(0), value);
237 :
238 0 : default:
239 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
240 : }
241 :
242 : // dummy
243 : return libMesh::RealVectorValue();
244 : }
245 : template <>
246 : RealVectorValue
247 0 : FE<2, MONOMIAL_VEC>::shape_deriv(const ElemType type,
248 : const Order order,
249 : const unsigned int i,
250 : const unsigned int j,
251 : const Point & p)
252 : {
253 0 : Real value = FE<2, MONOMIAL>::shape_deriv(type, order, i / 2, j, p);
254 :
255 0 : switch (i % 2)
256 : {
257 0 : case 0:
258 0 : return libMesh::RealVectorValue(value);
259 :
260 0 : case 1:
261 0 : return libMesh::RealVectorValue(Real(0), value);
262 :
263 0 : default:
264 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
265 : }
266 :
267 : // dummy
268 : return libMesh::RealVectorValue();
269 : }
270 :
271 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
272 : template <>
273 : RealVectorValue
274 0 : FE<2, MONOMIAL_VEC>::shape_second_deriv(const ElemType type,
275 : const Order order,
276 : const unsigned int i,
277 : const unsigned int j,
278 : const Point & p)
279 : {
280 0 : Real value = FE<2, MONOMIAL>::shape_second_deriv(type, order, i / 2, j, p);
281 :
282 0 : switch (i % 2)
283 : {
284 0 : case 0:
285 0 : return libMesh::RealVectorValue(value);
286 :
287 0 : case 1:
288 0 : return libMesh::RealVectorValue(Real(0), value);
289 :
290 0 : default:
291 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
292 : }
293 :
294 : // dummy
295 : return libMesh::RealVectorValue();
296 : }
297 :
298 : #endif
299 :
300 : // 3-D
301 : template <>
302 : RealVectorValue
303 0 : FE<3, MONOMIAL_VEC>::shape(const ElemType type,
304 : const Order order,
305 : const unsigned int i,
306 : const Point & p)
307 : {
308 0 : Real value = FE<3, MONOMIAL>::shape(type, order, i / 3, p);
309 :
310 0 : switch (i % 3)
311 : {
312 0 : case 0:
313 0 : return libMesh::RealVectorValue(value);
314 :
315 0 : case 1:
316 0 : return libMesh::RealVectorValue(Real(0), value);
317 :
318 0 : case 2:
319 0 : return libMesh::RealVectorValue(Real(0), Real(0), value);
320 :
321 0 : default:
322 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
323 : }
324 :
325 : // dummy
326 : return libMesh::RealVectorValue();
327 : }
328 : template <>
329 : RealVectorValue
330 0 : FE<3, MONOMIAL_VEC>::shape_deriv(const ElemType type,
331 : const Order order,
332 : const unsigned int i,
333 : const unsigned int j,
334 : const Point & p)
335 : {
336 0 : Real value = FE<3, MONOMIAL>::shape_deriv(type, order, i / 3, j, p);
337 :
338 0 : switch (i % 3)
339 : {
340 0 : case 0:
341 0 : return libMesh::RealVectorValue(value);
342 :
343 0 : case 1:
344 0 : return libMesh::RealVectorValue(Real(0), value);
345 :
346 0 : case 2:
347 0 : return libMesh::RealVectorValue(Real(0), Real(0), value);
348 :
349 0 : default:
350 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
351 : }
352 :
353 : // dummy
354 : return libMesh::RealVectorValue();
355 : }
356 :
357 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
358 :
359 : template <>
360 : RealVectorValue
361 0 : FE<3, MONOMIAL_VEC>::shape_second_deriv(const ElemType type,
362 : const Order order,
363 : const unsigned int i,
364 : const unsigned int j,
365 : const Point & p)
366 : {
367 0 : Real value = FE<3, MONOMIAL>::shape_second_deriv(type, order, i / 3, j, p);
368 :
369 0 : switch (i % 3)
370 : {
371 0 : case 0:
372 0 : return libMesh::RealVectorValue(value);
373 :
374 0 : case 1:
375 0 : return libMesh::RealVectorValue(Real(0), value);
376 :
377 0 : case 2:
378 0 : return libMesh::RealVectorValue(Real(0), Real(0), value);
379 :
380 0 : default:
381 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
382 : }
383 :
384 : // dummy
385 : return libMesh::RealVectorValue();
386 : }
387 :
388 : #endif
389 :
390 : // 0-D
391 : template <>
392 : RealVectorValue
393 0 : FE<0, MONOMIAL_VEC>::shape(const Elem * elem,
394 : const Order order,
395 : const unsigned int i,
396 : const Point & p,
397 : const bool add_p_level)
398 : {
399 : Real value =
400 0 : FE<0, MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i, p);
401 0 : return libMesh::RealVectorValue(value);
402 : }
403 : template <>
404 : RealVectorValue
405 0 : FE<0, MONOMIAL_VEC>::shape_deriv(const Elem * elem,
406 : const Order order,
407 : const unsigned int i,
408 : const unsigned int j,
409 : const Point & p,
410 : const bool add_p_level)
411 : {
412 0 : Real value = FE<0, MONOMIAL>::shape_deriv(
413 0 : elem->type(), order + add_p_level*elem->p_level(), i, j, p);
414 0 : return libMesh::RealVectorValue(value);
415 : }
416 :
417 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
418 :
419 : template <>
420 : RealVectorValue
421 0 : FE<0, MONOMIAL_VEC>::shape_second_deriv(const Elem * elem,
422 : const Order order,
423 : const unsigned int i,
424 : const unsigned int j,
425 : const Point & p,
426 : const bool add_p_level)
427 : {
428 0 : Real value = FE<0, MONOMIAL>::shape_second_deriv(
429 0 : elem->type(), order + add_p_level*elem->p_level(), i, j, p);
430 0 : return libMesh::RealVectorValue(value);
431 : }
432 :
433 : #endif
434 :
435 : // 1-D
436 : template <>
437 : RealVectorValue
438 0 : FE<1, MONOMIAL_VEC>::shape(const Elem * elem,
439 : const Order order,
440 : const unsigned int i,
441 : const Point & p,
442 : const bool add_p_level)
443 : {
444 : Real value =
445 0 : FE<1, MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i, p);
446 0 : return libMesh::RealVectorValue(value);
447 : }
448 : template <>
449 : RealVectorValue
450 0 : FE<1, MONOMIAL_VEC>::shape_deriv(const Elem * elem,
451 : const Order order,
452 : const unsigned int i,
453 : const unsigned int j,
454 : const Point & p,
455 : const bool add_p_level)
456 : {
457 0 : Real value = FE<1, MONOMIAL>::shape_deriv(
458 0 : elem->type(), order + add_p_level*elem->p_level(), i, j, p);
459 0 : return libMesh::RealVectorValue(value);
460 : }
461 :
462 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
463 : template <>
464 : RealVectorValue
465 0 : FE<1, MONOMIAL_VEC>::shape_second_deriv(const Elem * elem,
466 : const Order order,
467 : const unsigned int i,
468 : const unsigned int j,
469 : const Point & p,
470 : const bool add_p_level)
471 : {
472 0 : Real value = FE<1, MONOMIAL>::shape_second_deriv(
473 0 : elem->type(), order + add_p_level*elem->p_level(), i, j, p);
474 0 : return libMesh::RealVectorValue(value);
475 : }
476 :
477 : #endif
478 :
479 : // 2-D
480 : template <>
481 : RealVectorValue
482 1112430 : FE<2, MONOMIAL_VEC>::shape(const Elem * elem,
483 : const Order order,
484 : const unsigned int i,
485 : const Point & p,
486 : const bool add_p_level)
487 : {
488 : Real value =
489 1314690 : FE<2, MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i / 2, p);
490 :
491 1112430 : switch (i % 2)
492 : {
493 556215 : case 0:
494 50565 : return libMesh::RealVectorValue(value);
495 :
496 556215 : case 1:
497 50565 : return libMesh::RealVectorValue(Real(0), value);
498 :
499 0 : default:
500 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
501 : }
502 :
503 : // dummy
504 : return libMesh::RealVectorValue();
505 : }
506 : template <>
507 : RealVectorValue
508 1875920 : FE<2, MONOMIAL_VEC>::shape_deriv(const Elem * elem,
509 : const Order order,
510 : const unsigned int i,
511 : const unsigned int j,
512 : const Point & p,
513 : const bool add_p_level)
514 : {
515 1875920 : Real value = FE<2, MONOMIAL>::shape_deriv(
516 2045376 : elem->type(), order + add_p_level*elem->p_level(), i / 2, j, p);
517 :
518 1875920 : switch (i % 2)
519 : {
520 937960 : case 0:
521 84728 : return libMesh::RealVectorValue(value);
522 :
523 937960 : case 1:
524 84728 : return libMesh::RealVectorValue(Real(0), value);
525 :
526 0 : default:
527 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
528 : }
529 :
530 : // dummy
531 : return libMesh::RealVectorValue();
532 : }
533 :
534 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
535 : template <>
536 : RealVectorValue
537 0 : FE<2, MONOMIAL_VEC>::shape_second_deriv(const Elem * elem,
538 : const Order order,
539 : const unsigned int i,
540 : const unsigned int j,
541 : const Point & p,
542 : const bool add_p_level)
543 : {
544 0 : Real value = FE<2, MONOMIAL>::shape_second_deriv(
545 0 : elem->type(), order + add_p_level*elem->p_level(), i / 2, j, p);
546 :
547 0 : switch (i % 2)
548 : {
549 0 : case 0:
550 0 : return libMesh::RealVectorValue(value);
551 :
552 0 : case 1:
553 0 : return libMesh::RealVectorValue(Real(0), value);
554 :
555 0 : default:
556 0 : libmesh_error_msg("i%2 must be either 0 or 1!");
557 : }
558 :
559 : // dummy
560 : return libMesh::RealVectorValue();
561 : }
562 :
563 : #endif
564 :
565 : // 3-D
566 : template <>
567 : RealVectorValue
568 0 : FE<3, MONOMIAL_VEC>::shape(const Elem * elem,
569 : const Order order,
570 : const unsigned int i,
571 : const Point & p,
572 : const bool add_p_level)
573 : {
574 : Real value =
575 0 : FE<3, MONOMIAL>::shape(elem->type(), order + add_p_level*elem->p_level(), i / 3, p);
576 :
577 0 : switch (i % 3)
578 : {
579 0 : case 0:
580 0 : return libMesh::RealVectorValue(value);
581 :
582 0 : case 1:
583 0 : return libMesh::RealVectorValue(Real(0), value);
584 :
585 0 : case 2:
586 0 : return libMesh::RealVectorValue(Real(0), Real(0), value);
587 :
588 0 : default:
589 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
590 : }
591 :
592 : // dummy
593 : return libMesh::RealVectorValue();
594 : }
595 : template <>
596 : RealVectorValue
597 0 : FE<3, MONOMIAL_VEC>::shape_deriv(const Elem * elem,
598 : const Order order,
599 : const unsigned int i,
600 : const unsigned int j,
601 : const Point & p,
602 : const bool add_p_level)
603 : {
604 0 : Real value = FE<3, MONOMIAL>::shape_deriv(
605 0 : elem->type(), order + add_p_level*elem->p_level(), i / 3, j, p);
606 :
607 0 : switch (i % 3)
608 : {
609 0 : case 0:
610 0 : return libMesh::RealVectorValue(value);
611 :
612 0 : case 1:
613 0 : return libMesh::RealVectorValue(Real(0), value);
614 :
615 0 : case 2:
616 0 : return libMesh::RealVectorValue(Real(0), Real(0), value);
617 :
618 0 : default:
619 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
620 : }
621 :
622 : // dummy
623 : return libMesh::RealVectorValue();
624 : }
625 :
626 : #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
627 :
628 : template <>
629 : RealVectorValue
630 0 : FE<3, MONOMIAL_VEC>::shape_second_deriv(const Elem * elem,
631 : const Order order,
632 : const unsigned int i,
633 : const unsigned int j,
634 : const Point & p,
635 : const bool add_p_level)
636 : {
637 0 : Real value = FE<3, MONOMIAL>::shape_second_deriv(
638 0 : elem->type(), order + add_p_level*elem->p_level(), i / 3, j, p);
639 :
640 0 : switch (i % 3)
641 : {
642 0 : case 0:
643 0 : return libMesh::RealVectorValue(value);
644 :
645 0 : case 1:
646 0 : return libMesh::RealVectorValue(Real(0), value);
647 :
648 0 : case 2:
649 0 : return libMesh::RealVectorValue(Real(0), Real(0), value);
650 :
651 0 : default:
652 0 : libmesh_error_msg("i%3 must be 0, 1, or 2!");
653 : }
654 :
655 : // dummy
656 : return libMesh::RealVectorValue();
657 : }
658 :
659 : #endif
660 :
661 : // Do full-specialization for every dimension, instead
662 : // of explicit instantiation at the end of this function.
663 : // This could be macro-ified.
664 0 : template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs(const ElemType t, const Order o) { return FE<0, MONOMIAL>::n_dofs(t, o); }
665 0 : template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs(const ElemType t, const Order o) { return FE<1, MONOMIAL>::n_dofs(t, o); }
666 0 : template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs(const ElemType t, const Order o) { return 2 * FE<2, MONOMIAL>::n_dofs(t, o); }
667 0 : template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs(const ElemType t, const Order o) { return 3 * FE<3, MONOMIAL>::n_dofs(t, o); }
668 :
669 0 : template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs(const Elem * e, const Order o) { return FE<0, MONOMIAL>::n_dofs(e, o); }
670 0 : template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs(const Elem * e, const Order o) { return FE<1, MONOMIAL>::n_dofs(e, o); }
671 109815 : template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs(const Elem * e, const Order o) { return 2 * FE<2, MONOMIAL>::n_dofs(e, o); }
672 0 : template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs(const Elem * e, const Order o) { return 3 * FE<3, MONOMIAL>::n_dofs(e, o); }
673 :
674 0 : template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
675 0 : template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
676 343240 : template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
677 0 : template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int) { return 0; }
678 :
679 0 : template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
680 0 : template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
681 85292 : template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
682 0 : template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs_at_node(const Elem &, const Order, const unsigned int) { return 0; }
683 :
684 0 : template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<0, MONOMIAL>::n_dofs_per_elem(t, o); }
685 0 : template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<1, MONOMIAL>::n_dofs_per_elem(t, o); }
686 0 : template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 2 * FE<2, MONOMIAL>::n_dofs_per_elem(t, o); }
687 0 : template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 3 * FE<3, MONOMIAL>::n_dofs_per_elem(t, o); }
688 :
689 0 : template <> unsigned int FE<0, MONOMIAL_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<0, MONOMIAL>::n_dofs_per_elem(e, o); }
690 0 : template <> unsigned int FE<1, MONOMIAL_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<1, MONOMIAL>::n_dofs_per_elem(e, o); }
691 96494 : template <> unsigned int FE<2, MONOMIAL_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 2 * FE<2, MONOMIAL>::n_dofs_per_elem(e, o); }
692 0 : template <> unsigned int FE<3, MONOMIAL_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 3 * FE<3, MONOMIAL>::n_dofs_per_elem(e, o); }
693 :
694 : // Monomial FEMs are always C^0 continuous
695 : template <>
696 : FEContinuity
697 0 : FE<0, MONOMIAL_VEC>::get_continuity() const
698 : {
699 0 : return DISCONTINUOUS;
700 : }
701 : template <>
702 : FEContinuity
703 0 : FE<1, MONOMIAL_VEC>::get_continuity() const
704 : {
705 0 : return DISCONTINUOUS;
706 : }
707 : template <>
708 : FEContinuity
709 355 : FE<2, MONOMIAL_VEC>::get_continuity() const
710 : {
711 355 : return DISCONTINUOUS;
712 : }
713 : template <>
714 : FEContinuity
715 0 : FE<3, MONOMIAL_VEC>::get_continuity() const
716 : {
717 0 : return DISCONTINUOUS;
718 : }
719 :
720 : // Monomial FEMs are not hierarchic
721 : template <>
722 : bool
723 0 : FE<0, MONOMIAL_VEC>::is_hierarchic() const
724 : {
725 0 : return true;
726 : }
727 : template <>
728 : bool
729 0 : FE<1, MONOMIAL_VEC>::is_hierarchic() const
730 : {
731 0 : return true;
732 : }
733 : template <>
734 : bool
735 0 : FE<2, MONOMIAL_VEC>::is_hierarchic() const
736 : {
737 0 : return true;
738 : }
739 : template <>
740 : bool
741 0 : FE<3, MONOMIAL_VEC>::is_hierarchic() const
742 : {
743 0 : return true;
744 : }
745 :
746 : // Monomial FEM shapes do not need reinit (is this always true?)
747 : template <>
748 : bool
749 0 : FE<0, MONOMIAL_VEC>::shapes_need_reinit() const
750 : {
751 0 : return false;
752 : }
753 : template <>
754 : bool
755 0 : FE<1, MONOMIAL_VEC>::shapes_need_reinit() const
756 : {
757 0 : return false;
758 : }
759 : template <>
760 : bool
761 102970 : FE<2, MONOMIAL_VEC>::shapes_need_reinit() const
762 : {
763 102970 : return false;
764 : }
765 : template <>
766 : bool
767 0 : FE<3, MONOMIAL_VEC>::shapes_need_reinit() const
768 : {
769 0 : return false;
770 : }
771 :
772 : // we only need instantiations of this function for
773 : // Dim==2 and 3. There are no constraints for discontinuous elements, so
774 : // we do nothing.
775 : #ifdef LIBMESH_ENABLE_AMR
776 : template <>
777 : void
778 0 : FE<2, MONOMIAL_VEC>::compute_constraints(DofConstraints &,
779 : DofMap &,
780 : const unsigned int,
781 : const Elem *)
782 : {
783 0 : }
784 :
785 : template <>
786 : void
787 0 : FE<3, MONOMIAL_VEC>::compute_constraints(DofConstraints &,
788 : DofMap &,
789 : const unsigned int,
790 : const Elem *)
791 : {
792 0 : }
793 : #endif // LIBMESH_ENABLE_AMR
794 :
795 : } // namespace libMesh
|