libMesh
fe_monomial_vec.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2019 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/fe.h"
21 #include "libmesh/fe_interface.h"
22 #include "libmesh/elem.h"
23 #include "libmesh/tensor_value.h"
24 
25 namespace libMesh
26 {
27 
28 // ------------------------------------------------------------
29 // Vector monomial specific implementations
30 
31 // Anonymous namespace for local helper functions
32 namespace
33 {
34 void
35 monomial_vec_nodal_soln(const Elem * elem,
36  const Order order,
37  const std::vector<Number> & elem_soln,
38  const int dim,
39  std::vector<Number> & nodal_soln)
40 {
41  const unsigned int n_nodes = elem->n_nodes();
42 
43  const ElemType elem_type = elem->type();
44 
45  nodal_soln.resize(dim * n_nodes);
46 
47  const Order totalorder = static_cast<Order>(order + elem->p_level());
48 
49  switch (totalorder)
50  {
51  // Constant shape functions
52  case CONSTANT:
53  {
54  switch (dim)
55  {
56  case 2:
57  {
58  libmesh_assert_equal_to(elem_soln.size(), 2);
59  for (unsigned int n = 0; n < n_nodes; n++)
60  {
61  nodal_soln[2 * n] = elem_soln[0];
62  nodal_soln[1 + 2 * n] = elem_soln[1];
63  }
64  return;
65  }
66  case 3:
67  {
68  libmesh_assert_equal_to(elem_soln.size(), 3);
69  for (unsigned int n = 0; n < n_nodes; n++)
70  {
71  nodal_soln[3 * n] = elem_soln[0];
72  nodal_soln[1 + 3 * n] = elem_soln[1];
73  nodal_soln[2 + 3 * n] = elem_soln[2];
74  }
75  return;
76  }
77  default:
78  libmesh_error_msg(
79  "The monomial_vec_nodal_soln helper should only be called for 2 and 3 dimensions");
80  }
81  }
82 
83  // For other orders, do interpolation at the nodes
84  // explicitly.
85  default:
86  {
87  // FEType object to be passed to various FEInterface functions below.
88  FEType fe_type(totalorder, MONOMIAL);
89 
90  const unsigned int n_sf = FEInterface::n_shape_functions(dim, fe_type, elem_type);
91 
92  std::vector<Point> refspace_nodes;
93  FEBase::get_refspace_nodes(elem_type, refspace_nodes);
94  libmesh_assert_equal_to(refspace_nodes.size(), n_nodes);
95  libmesh_assert_equal_to(elem_soln.size(), n_sf * dim);
96 
97  for (unsigned int d = 0; d < static_cast<unsigned int>(dim); d++)
98  for (unsigned int n = 0; n < n_nodes; n++)
99  {
100 
101  // Zero before summation
102  nodal_soln[d + dim * n] = 0;
103 
104  // u_i = Sum (alpha_i phi_i)
105  for (unsigned int i = 0; i < n_sf; i++)
106  nodal_soln[d + dim * n] += elem_soln[d + dim * i] *
107  FEInterface::shape(dim, fe_type, elem, i, refspace_nodes[n]);
108  }
109 
110  return;
111  } // default
112  } // switch
113 }
114 } // anonymous namespace
115 
116 // Do full-specialization for every dimension, instead
117 // of explicit instantiation at the end of this file.
118 // This could be macro-ified so that it fits on one line...
119 template <>
120 void
122  const Order order,
123  const std::vector<Number> & elem_soln,
124  std::vector<Number> & nodal_soln)
125 {
126  FE<0, MONOMIAL>::nodal_soln(elem, order, elem_soln, nodal_soln);
127 }
128 
129 template <>
130 void
132  const Order order,
133  const std::vector<Number> & elem_soln,
134  std::vector<Number> & nodal_soln)
135 {
136  FE<1, MONOMIAL>::nodal_soln(elem, order, elem_soln, nodal_soln);
137 }
138 
139 template <>
140 void
142  const Order order,
143  const std::vector<Number> & elem_soln,
144  std::vector<Number> & nodal_soln)
145 {
146  monomial_vec_nodal_soln(elem, order, elem_soln, 2 /*dimension*/, nodal_soln);
147 }
148 
149 template <>
150 void
152  const Order order,
153  const std::vector<Number> & elem_soln,
154  std::vector<Number> & nodal_soln)
155 {
156  monomial_vec_nodal_soln(elem, order, elem_soln, 3 /*dimension*/, nodal_soln);
157 }
158 
159 // Specialize for shape function routines by leveraging scalar MONOMIAL elements
160 
161 // 0-D
162 template <>
165  const Order order,
166  const unsigned int i,
167  const Point & p)
168 {
169  Real value = FE<0, MONOMIAL>::shape(type, order, i, p);
171 }
172 template <>
175  const Order order,
176  const unsigned int i,
177  const unsigned int j,
178  const Point & p)
179 {
180  Real value = FE<0, MONOMIAL>::shape_deriv(type, order, i, j, p);
182 }
183 
184 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
185 
186 template <>
189  const Order order,
190  const unsigned int i,
191  const unsigned int j,
192  const Point & p)
193 {
194  Real value = FE<0, MONOMIAL>::shape_second_deriv(type, order, i, j, p);
196 }
197 #endif
198 
199 // 1-D
200 template <>
203  const Order order,
204  const unsigned int i,
205  const Point & p)
206 {
207  Real value = FE<1, MONOMIAL>::shape(type, order, i, p);
209 }
210 template <>
213  const Order order,
214  const unsigned int i,
215  const unsigned int j,
216  const Point & p)
217 {
218  Real value = FE<1, MONOMIAL>::shape_deriv(type, order, i, j, p);
220 }
221 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
222 template <>
225  const Order order,
226  const unsigned int i,
227  const unsigned int j,
228  const Point & p)
229 {
230  Real value = FE<1, MONOMIAL>::shape_second_deriv(type, order, i, j, p);
232 }
233 
234 #endif
235 
236 // 2-D
237 template <>
240  const Order order,
241  const unsigned int i,
242  const Point & p)
243 {
244  Real value = FE<2, MONOMIAL>::shape(type, order, i / 2, p);
245 
246  switch (i % 2)
247  {
248  case 0:
250 
251  case 1:
252  return libMesh::RealVectorValue(Real(0), value);
253 
254  default:
255  libmesh_error_msg("i%2 must be either 0 or 1!");
256  }
257 
258  // dummy
259  return libMesh::RealVectorValue();
260 }
261 template <>
264  const Order order,
265  const unsigned int i,
266  const unsigned int j,
267  const Point & p)
268 {
269  Real value = FE<2, MONOMIAL>::shape_deriv(type, order, i / 2, j, p);
270 
271  switch (i % 2)
272  {
273  case 0:
275 
276  case 1:
277  return libMesh::RealVectorValue(Real(0), value);
278 
279  default:
280  libmesh_error_msg("i%2 must be either 0 or 1!");
281  }
282 
283  // dummy
284  return libMesh::RealVectorValue();
285 }
286 
287 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
288 template <>
291  const Order order,
292  const unsigned int i,
293  const unsigned int j,
294  const Point & p)
295 {
296  Real value = FE<2, MONOMIAL>::shape_second_deriv(type, order, i / 2, j, p);
297 
298  switch (i % 2)
299  {
300  case 0:
302 
303  case 1:
304  return libMesh::RealVectorValue(Real(0), value);
305 
306  default:
307  libmesh_error_msg("i%2 must be either 0 or 1!");
308  }
309 
310  // dummy
311  return libMesh::RealVectorValue();
312 }
313 
314 #endif
315 
316 // 3-D
317 template <>
320  const Order order,
321  const unsigned int i,
322  const Point & p)
323 {
324  Real value = FE<3, MONOMIAL>::shape(type, order, i / 3, p);
325 
326  switch (i % 3)
327  {
328  case 0:
330 
331  case 1:
332  return libMesh::RealVectorValue(Real(0), value);
333 
334  case 2:
335  return libMesh::RealVectorValue(Real(0), Real(0), value);
336 
337  default:
338  libmesh_error_msg("i%3 must be 0, 1, or 2!");
339  }
340 
341  // dummy
342  return libMesh::RealVectorValue();
343 }
344 template <>
347  const Order order,
348  const unsigned int i,
349  const unsigned int j,
350  const Point & p)
351 {
352  Real value = FE<3, MONOMIAL>::shape_deriv(type, order, i / 3, j, p);
353 
354  switch (i % 3)
355  {
356  case 0:
358 
359  case 1:
360  return libMesh::RealVectorValue(Real(0), value);
361 
362  case 2:
363  return libMesh::RealVectorValue(Real(0), Real(0), value);
364 
365  default:
366  libmesh_error_msg("i%3 must be 0, 1, or 2!");
367  }
368 
369  // dummy
370  return libMesh::RealVectorValue();
371 }
372 
373 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
374 
375 template <>
378  const Order order,
379  const unsigned int i,
380  const unsigned int j,
381  const Point & p)
382 {
383  Real value = FE<3, MONOMIAL>::shape_second_deriv(type, order, i / 3, j, p);
384 
385  switch (i % 3)
386  {
387  case 0:
389 
390  case 1:
391  return libMesh::RealVectorValue(Real(0), value);
392 
393  case 2:
394  return libMesh::RealVectorValue(Real(0), Real(0), value);
395 
396  default:
397  libmesh_error_msg("i%3 must be 0, 1, or 2!");
398  }
399 
400  // dummy
401  return libMesh::RealVectorValue();
402 }
403 
404 #endif
405 
406 // 0-D
407 template <>
410  const Order order,
411  const unsigned int i,
412  const Point & p,
413  const bool add_p_level)
414 {
415  Real value =
416  FE<0, MONOMIAL>::shape(elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, p);
418 }
419 template <>
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 {
429  elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, j, p);
431 }
432 
433 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
434 
435 template <>
438  const Order order,
439  const unsigned int i,
440  const unsigned int j,
441  const Point & p,
442  const bool add_p_level)
443 {
445  elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, j, p);
447 }
448 
449 #endif
450 
451 // 1-D
452 template <>
455  const Order order,
456  const unsigned int i,
457  const Point & p,
458  const bool add_p_level)
459 {
460  Real value =
461  FE<1, MONOMIAL>::shape(elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, p);
463 }
464 template <>
467  const Order order,
468  const unsigned int i,
469  const unsigned int j,
470  const Point & p,
471  const bool add_p_level)
472 {
474  elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, j, p);
476 }
477 
478 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
479 template <>
482  const Order order,
483  const unsigned int i,
484  const unsigned int j,
485  const Point & p,
486  const bool add_p_level)
487 {
489  elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i, j, p);
491 }
492 
493 #endif
494 
495 // 2-D
496 template <>
499  const Order order,
500  const unsigned int i,
501  const Point & p,
502  const bool add_p_level)
503 {
504  Real value =
505  FE<2, MONOMIAL>::shape(elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i / 2, p);
506 
507  switch (i % 2)
508  {
509  case 0:
511 
512  case 1:
513  return libMesh::RealVectorValue(Real(0), value);
514 
515  default:
516  libmesh_error_msg("i%2 must be either 0 or 1!");
517  }
518 
519  // dummy
520  return libMesh::RealVectorValue();
521 }
522 template <>
525  const Order order,
526  const unsigned int i,
527  const unsigned int j,
528  const Point & p,
529  const bool add_p_level)
530 {
532  elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i / 2, j, p);
533 
534  switch (i % 2)
535  {
536  case 0:
538 
539  case 1:
540  return libMesh::RealVectorValue(Real(0), value);
541 
542  default:
543  libmesh_error_msg("i%2 must be either 0 or 1!");
544  }
545 
546  // dummy
547  return libMesh::RealVectorValue();
548 }
549 
550 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
551 template <>
554  const Order order,
555  const unsigned int i,
556  const unsigned int j,
557  const Point & p,
558  const bool add_p_level)
559 {
561  elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i / 2, j, p);
562 
563  switch (i % 2)
564  {
565  case 0:
567 
568  case 1:
569  return libMesh::RealVectorValue(Real(0), value);
570 
571  default:
572  libmesh_error_msg("i%2 must be either 0 or 1!");
573  }
574 
575  // dummy
576  return libMesh::RealVectorValue();
577 }
578 
579 #endif
580 
581 // 3-D
582 template <>
585  const Order order,
586  const unsigned int i,
587  const Point & p,
588  const bool add_p_level)
589 {
590  Real value =
591  FE<3, MONOMIAL>::shape(elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i / 3, p);
592 
593  switch (i % 3)
594  {
595  case 0:
597 
598  case 1:
599  return libMesh::RealVectorValue(Real(0), value);
600 
601  case 2:
602  return libMesh::RealVectorValue(Real(0), Real(0), value);
603 
604  default:
605  libmesh_error_msg("i%3 must be 0, 1, or 2!");
606  }
607 
608  // dummy
609  return libMesh::RealVectorValue();
610 }
611 template <>
614  const Order order,
615  const unsigned int i,
616  const unsigned int j,
617  const Point & p,
618  const bool add_p_level)
619 {
621  elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i / 3, j, p);
622 
623  switch (i % 3)
624  {
625  case 0:
627 
628  case 1:
629  return libMesh::RealVectorValue(Real(0), value);
630 
631  case 2:
632  return libMesh::RealVectorValue(Real(0), Real(0), value);
633 
634  default:
635  libmesh_error_msg("i%3 must be 0, 1, or 2!");
636  }
637 
638  // dummy
639  return libMesh::RealVectorValue();
640 }
641 
642 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES
643 
644 template <>
647  const Order order,
648  const unsigned int i,
649  const unsigned int j,
650  const Point & p,
651  const bool add_p_level)
652 {
654  elem->type(), static_cast<Order>(order + add_p_level * elem->p_level()), i / 3, j, p);
655 
656  switch (i % 3)
657  {
658  case 0:
660 
661  case 1:
662  return libMesh::RealVectorValue(Real(0), value);
663 
664  case 2:
665  return libMesh::RealVectorValue(Real(0), Real(0), value);
666 
667  default:
668  libmesh_error_msg("i%3 must be 0, 1, or 2!");
669  }
670 
671  // dummy
672  return libMesh::RealVectorValue();
673 }
674 
675 #endif
676 
677 // Do full-specialization for every dimension, instead
678 // of explicit instantiation at the end of this function.
679 // This could be macro-ified.
680 template <>
681 unsigned int
683 {
684  return FE<0, MONOMIAL>::n_dofs(t, o);
685 }
686 template <>
687 unsigned int
689 {
690  return FE<1, MONOMIAL>::n_dofs(t, o);
691 }
692 template <>
693 unsigned int
695 {
696  return 2 * FE<2, MONOMIAL>::n_dofs(t, o);
697 }
698 template <>
699 unsigned int
701 {
702  return 3 * FE<3, MONOMIAL>::n_dofs(t, o);
703 }
704 
705 template <>
706 unsigned int
707 FE<0, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int)
708 {
709  return 0;
710 }
711 template <>
712 unsigned int
713 FE<1, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int)
714 {
715  return 0;
716 }
717 template <>
718 unsigned int
719 FE<2, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int)
720 {
721  return 0;
722 }
723 template <>
724 unsigned int
725 FE<3, MONOMIAL_VEC>::n_dofs_at_node(const ElemType, const Order, const unsigned int)
726 {
727  return 0;
728 }
729 
730 template <>
731 unsigned int
733 {
735 }
736 template <>
737 unsigned int
739 {
741 }
742 template <>
743 unsigned int
745 {
746  return 2 * FE<2, MONOMIAL>::n_dofs_per_elem(t, o);
747 }
748 template <>
749 unsigned int
751 {
752  return 3 * FE<3, MONOMIAL>::n_dofs_per_elem(t, o);
753 }
754 
755 // Monomial FEMs are always C^0 continuous
756 template <>
759 {
760  return DISCONTINUOUS;
761 }
762 template <>
765 {
766  return DISCONTINUOUS;
767 }
768 template <>
771 {
772  return DISCONTINUOUS;
773 }
774 template <>
777 {
778  return DISCONTINUOUS;
779 }
780 
781 // Monomial FEMs are not hierarchic
782 template <>
783 bool
785 {
786  return true;
787 }
788 template <>
789 bool
791 {
792  return true;
793 }
794 template <>
795 bool
797 {
798  return true;
799 }
800 template <>
801 bool
803 {
804  return true;
805 }
806 
807 // Monomial FEM shapes do not need reinit (is this always true?)
808 template <>
809 bool
811 {
812  return false;
813 }
814 template <>
815 bool
817 {
818  return false;
819 }
820 template <>
821 bool
823 {
824  return false;
825 }
826 template <>
827 bool
829 {
830  return false;
831 }
832 
833 // we only need instantiations of this function for
834 // Dim==2 and 3. There are no constraints for discontinuous elements, so
835 // we do nothing.
836 #ifdef LIBMESH_ENABLE_AMR
837 template <>
838 void
840  DofMap &,
841  const unsigned int,
842  const Elem *)
843 {
844 }
845 
846 template <>
847 void
849  DofMap &,
850  const unsigned int,
851  const Elem *)
852 {
853 }
854 #endif // LIBMESH_ENABLE_AMR
855 
856 } // namespace libMesh
libMesh::DofConstraints
The constraint matrix storage format.
Definition: dof_map.h:105
libMesh::FE::shape_second_deriv
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::FEInterface::n_shape_functions
static unsigned int n_shape_functions(const unsigned int dim, const FEType &fe_t, const ElemType t)
Definition: fe_interface.C:446
libMesh
The libMesh namespace provides an interface to certain functionality in the library.
Definition: factoryfunction.C:55
libMesh::Order
Order
Definition: enum_order.h:40
libMesh::Elem::p_level
unsigned int p_level() const
Definition: elem.h:2512
libMesh::RealVectorValue
VectorValue< Real > RealVectorValue
Useful typedefs to allow transparent switching between Real and Complex data types.
Definition: hp_coarsentest.h:46
dim
unsigned int dim
Definition: adaptivity_ex3.C:113
libMesh::VectorValue
This class defines a vector in LIBMESH_DIM dimensional Real or Complex space.
Definition: exact_solution.h:54
libMesh::FE::shape
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
libMesh::FE::shape_deriv
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
libMesh::Point
A Point defines a location in LIBMESH_DIM dimensional Real space.
Definition: point.h:38
libMesh::DISCONTINUOUS
Definition: enum_fe_family.h:78
libMesh::FE::compute_constraints
static void compute_constraints(DofConstraints &constraints, DofMap &dof_map, const unsigned int variable_number, const Elem *elem)
Computes the constraint matrix contributions (for non-conforming adapted meshes) corresponding to var...
libMesh::CONSTANT
Definition: enum_order.h:41
n_nodes
const dof_id_type n_nodes
Definition: tecplot_io.C:68
libMesh::MONOMIAL
Definition: enum_fe_family.h:39
libMesh::FE::shapes_need_reinit
virtual bool shapes_need_reinit() const override
libMesh::FE::n_dofs_at_node
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
libMesh::FE::nodal_soln
static void nodal_soln(const Elem *elem, const Order o, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln)
Build the nodal soln from the element soln.
libMesh::FE::get_continuity
virtual FEContinuity get_continuity() const override
value
static const bool value
Definition: xdr_io.C:56
libMesh::DofMap
This class handles the numbering of degrees of freedom on a mesh.
Definition: dof_map.h:176
libMesh::FE::n_dofs_per_elem
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
libMesh::FE::is_hierarchic
virtual bool is_hierarchic() const override
libMesh::Elem
This is the base class from which all geometric element types are derived.
Definition: elem.h:100
libMesh::FEAbstract::get_refspace_nodes
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
Definition: fe_abstract.C:308
libMesh::FE::n_dofs
static unsigned int n_dofs(const ElemType t, const Order o)
libMesh::Real
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
Definition: libmesh_common.h:121
libMesh::FEContinuity
FEContinuity
Definition: enum_fe_family.h:77
libMesh::Elem::type
virtual ElemType type() const =0
libMesh::ElemType
ElemType
Defines an enum for geometric element types.
Definition: enum_elem_type.h:33
libMesh::FEInterface::shape
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
Definition: fe_interface.C:687