19 #include "libmesh/fe.h" 20 #include "libmesh/libmesh_logging.h" 21 #include "libmesh/enum_elem_type.h" 22 #include "libmesh/boundary_info.h" 23 #include "libmesh/mesh_base.h" 24 #include "libmesh/dense_matrix.h" 25 #include "libmesh/dense_vector.h" 26 #include "libmesh/dof_map.h" 27 #include "libmesh/elem.h" 28 #include "libmesh/fe_interface.h" 29 #include "libmesh/numeric_vector.h" 30 #include "libmesh/periodic_boundaries.h" 31 #include "libmesh/periodic_boundary.h" 32 #include "libmesh/quadrature.h" 33 #include "libmesh/quadrature_gauss.h" 34 #include "libmesh/remote_elem.h" 35 #include "libmesh/tensor_value.h" 36 #include "libmesh/threads.h" 37 #include "libmesh/enum_to_string.h" 39 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 40 #include "libmesh/inf_fe.h" 41 #include "libmesh/fe_interface_macros.h" 49 _fe_map(
FEMap::build(fet) ),
51 calculations_started(false),
52 calculate_dual(false),
53 calculate_default_dual_coeff(true),
54 calculate_nothing(false),
57 calculate_dphi(false),
58 calculate_d2phi(false),
59 calculate_curl_phi(false),
60 calculate_div_phi(false),
61 calculate_dphiref(false),
68 shapes_on_quadrature(false),
70 _add_p_level_in_reinit(true)
89 return std::make_unique<FE<0,CLOUGH>>(fet);
92 return std::make_unique<FE<0,HERMITE>>(fet);
95 return std::make_unique<FE<0,LAGRANGE>>(fet);
98 return std::make_unique<FE<0,LAGRANGE_VEC>>(fet);
101 return std::make_unique<FE<0,L2_LAGRANGE>>(fet);
104 return std::make_unique<FE<0,L2_LAGRANGE_VEC>>(fet);
107 return std::make_unique<FE<0,HIERARCHIC_VEC>>(fet);
110 return std::make_unique<FE<0,HIERARCHIC>>(fet);
113 return std::make_unique<FE<0,L2_HIERARCHIC>>(fet);
116 return std::make_unique<FE<0,L2_HIERARCHIC_VEC>>(fet);
119 return std::make_unique<FE<0,SIDE_HIERARCHIC>>(fet);
122 return std::make_unique<FE<0,MONOMIAL>>(fet);
125 return std::make_unique<FE<0,MONOMIAL_VEC>>(fet);
127 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 129 return std::make_unique<FE<0,SZABAB>>(fet);
132 return std::make_unique<FE<0,BERNSTEIN>>(fet);
135 return std::make_unique<FE<0,RATIONAL_BERNSTEIN>>(fet);
139 return std::make_unique<FEXYZ<0>>(fet);
142 return std::make_unique<FEScalar<0>>(fet);
145 return std::make_unique<FENedelecOne<0>>(fet);
148 return std::make_unique<FERaviartThomas<0>>(fet);
151 return std::make_unique<FEL2RaviartThomas<0>>(fet);
154 return std::make_unique<FESubdivision>(fet);
166 return std::make_unique<FE<1,CLOUGH>>(fet);
169 return std::make_unique<FE<1,HERMITE>>(fet);
172 return std::make_unique<FE<1,LAGRANGE>>(fet);
175 return std::make_unique<FE<1,LAGRANGE_VEC>>(fet);
178 return std::make_unique<FE<1,L2_LAGRANGE>>(fet);
181 return std::make_unique<FE<1,L2_LAGRANGE_VEC>>(fet);
184 return std::make_unique<FE<1,HIERARCHIC_VEC>>(fet);
187 return std::make_unique<FE<1,HIERARCHIC>>(fet);
190 return std::make_unique<FE<1,L2_HIERARCHIC>>(fet);
193 return std::make_unique<FE<1,L2_HIERARCHIC_VEC>>(fet);
196 return std::make_unique<FE<1,SIDE_HIERARCHIC>>(fet);
199 return std::make_unique<FE<1,MONOMIAL>>(fet);
202 return std::make_unique<FE<1,MONOMIAL_VEC>>(fet);
204 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 206 return std::make_unique<FE<1,SZABAB>>(fet);
209 return std::make_unique<FE<1,BERNSTEIN>>(fet);
212 return std::make_unique<FE<1,RATIONAL_BERNSTEIN>>(fet);
216 return std::make_unique<FEXYZ<1>>(fet);
219 return std::make_unique<FEScalar<1>>(fet);
222 return std::make_unique<FENedelecOne<1>>(fet);
225 return std::make_unique<FERaviartThomas<1>>(fet);
228 return std::make_unique<FEL2RaviartThomas<1>>(fet);
231 return std::make_unique<FESubdivision>(fet);
245 return std::make_unique<FE<2,CLOUGH>>(fet);
248 return std::make_unique<FE<2,HERMITE>>(fet);
251 return std::make_unique<FE<2,LAGRANGE>>(fet);
254 return std::make_unique<FE<2,LAGRANGE_VEC>>(fet);
257 return std::make_unique<FE<2,L2_LAGRANGE>>(fet);
260 return std::make_unique<FE<2,L2_LAGRANGE_VEC>>(fet);
263 return std::make_unique<FE<2,HIERARCHIC_VEC>>(fet);
266 return std::make_unique<FE<2,HIERARCHIC>>(fet);
269 return std::make_unique<FE<2,L2_HIERARCHIC>>(fet);
272 return std::make_unique<FE<2,L2_HIERARCHIC_VEC>>(fet);
275 return std::make_unique<FE<2,SIDE_HIERARCHIC>>(fet);
278 return std::make_unique<FE<2,MONOMIAL>>(fet);
281 return std::make_unique<FE<2,MONOMIAL_VEC>>(fet);
283 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 285 return std::make_unique<FE<2,SZABAB>>(fet);
288 return std::make_unique<FE<2,BERNSTEIN>>(fet);
291 return std::make_unique<FE<2,RATIONAL_BERNSTEIN>>(fet);
295 return std::make_unique<FEXYZ<2>>(fet);
298 return std::make_unique<FEScalar<2>>(fet);
301 return std::make_unique<FENedelecOne<2>>(fet);
304 return std::make_unique<FERaviartThomas<2>>(fet);
307 return std::make_unique<FEL2RaviartThomas<2>>(fet);
310 return std::make_unique<FESubdivision>(fet);
324 libmesh_error_msg(
"ERROR: Clough-Tocher elements currently only support 1D and 2D");
327 return std::make_unique<FE<3,HERMITE>>(fet);
330 return std::make_unique<FE<3,LAGRANGE>>(fet);
333 return std::make_unique<FE<3,LAGRANGE_VEC>>(fet);
336 return std::make_unique<FE<3,L2_LAGRANGE>>(fet);
339 return std::make_unique<FE<3,L2_LAGRANGE_VEC>>(fet);
342 return std::make_unique<FE<3,HIERARCHIC_VEC>>(fet);
345 return std::make_unique<FE<3,HIERARCHIC>>(fet);
348 return std::make_unique<FE<3,L2_HIERARCHIC>>(fet);
351 return std::make_unique<FE<3,L2_HIERARCHIC_VEC>>(fet);
354 return std::make_unique<FE<3,SIDE_HIERARCHIC>>(fet);
357 return std::make_unique<FE<3,MONOMIAL>>(fet);
360 return std::make_unique<FE<3,MONOMIAL_VEC>>(fet);
362 #ifdef LIBMESH_ENABLE_HIGHER_ORDER_SHAPES 364 return std::make_unique<FE<3,SZABAB>>(fet);
367 return std::make_unique<FE<3,BERNSTEIN>>(fet);
370 return std::make_unique<FE<3,RATIONAL_BERNSTEIN>>(fet);
374 return std::make_unique<FEXYZ<3>>(fet);
377 return std::make_unique<FEScalar<3>>(fet);
380 return std::make_unique<FENedelecOne<3>>(fet);
383 return std::make_unique<FERaviartThomas<3>>(fet);
386 return std::make_unique<FEL2RaviartThomas<3>>(fet);
394 libmesh_error_msg(
"Invalid dimension dim = " <<
dim);
404 libmesh_error_msg(
"Number of nodes is not well-defined for " <<
412 nodes[0] =
Point (0.,0.,0.);
417 nodes[2] =
Point (0.,0.,0.);
418 libmesh_fallthrough();
422 nodes[0] =
Point (-1.,0.,0.);
423 nodes[1] =
Point (1.,0.,0.);
428 nodes[0] =
Point (-1.,0.,0.);
429 nodes[1] =
Point (1.,0.,0.);
430 nodes[2] =
Point (-1./3.,0.,0.);
431 nodes[3] -
Point (1./3.,0.,0.);
436 nodes[6] =
Point (1./3.,1./3.,0.);
437 libmesh_fallthrough();
441 nodes[3] =
Point (.5,0.,0.);
442 nodes[4] =
Point (.5,.5,0.);
443 nodes[5] =
Point (0.,.5,0.);
444 libmesh_fallthrough();
449 nodes[0] =
Point (0.,0.,0.);
450 nodes[1] =
Point (1.,0.,0.);
451 nodes[2] =
Point (0.,1.,0.);
457 nodes[8] =
Point (0.,0.,0.);
458 libmesh_fallthrough();
463 nodes[4] =
Point (0.,-1.,0.);
464 nodes[5] =
Point (1.,0.,0.);
465 nodes[6] =
Point (0.,1.,0.);
466 nodes[7] =
Point (-1.,0.,0.);
467 libmesh_fallthrough();
472 nodes[0] =
Point (-1.,-1.,0.);
473 nodes[1] =
Point (1.,-1.,0.);
474 nodes[2] =
Point (1.,1.,0.);
475 nodes[3] =
Point (-1.,1.,0.);
484 libmesh_fallthrough();
488 nodes[4] =
Point (.5,0.,0.);
489 nodes[5] =
Point (.5,.5,0.);
490 nodes[6] =
Point (0.,.5,0.);
491 nodes[7] =
Point (0.,0.,.5);
492 nodes[8] =
Point (.5,0.,.5);
493 nodes[9] =
Point (0.,.5,.5);
494 libmesh_fallthrough();
498 nodes[0] =
Point (0.,0.,0.);
499 nodes[1] =
Point (1.,0.,0.);
500 nodes[2] =
Point (0.,1.,0.);
501 nodes[3] =
Point (0.,0.,1.);
506 nodes[20] =
Point (0.,0.,-1.);
507 nodes[21] =
Point (0.,-1.,0.);
508 nodes[22] =
Point (1.,0.,0.);
509 nodes[23] =
Point (0.,1.,0.);
510 nodes[24] =
Point (-1.,0.,0.);
511 nodes[25] =
Point (0.,0.,1.);
512 nodes[26] =
Point (0.,0.,0.);
513 libmesh_fallthrough();
517 nodes[8] =
Point (0.,-1.,-1.);
518 nodes[9] =
Point (1.,0.,-1.);
519 nodes[10] =
Point (0.,1.,-1.);
520 nodes[11] =
Point (-1.,0.,-1.);
521 nodes[12] =
Point (-1.,-1.,0.);
522 nodes[13] =
Point (1.,-1.,0.);
523 nodes[14] =
Point (1.,1.,0.);
524 nodes[15] =
Point (-1.,1.,0.);
525 nodes[16] =
Point (0.,-1.,1.);
526 nodes[17] =
Point (1.,0.,1.);
527 nodes[18] =
Point (0.,1.,1.);
528 nodes[19] =
Point (-1.,0.,1.);
529 libmesh_fallthrough();
533 nodes[0] =
Point (-1.,-1.,-1.);
534 nodes[1] =
Point (1.,-1.,-1.);
535 nodes[2] =
Point (1.,1.,-1.);
536 nodes[3] =
Point (-1.,1.,-1.);
537 nodes[4] =
Point (-1.,-1.,1.);
538 nodes[5] =
Point (1.,-1.,1.);
539 nodes[6] =
Point (1.,1.,1.);
540 nodes[7] =
Point (-1.,1.,1.);
546 libmesh_fallthrough();
552 libmesh_fallthrough();
556 nodes[15] =
Point (.5,0.,0.);
557 nodes[16] =
Point (.5,.5,0.);
558 nodes[17] =
Point (0.,.5,0.);
559 libmesh_fallthrough();
563 nodes[6] =
Point (.5,0.,-1.);
564 nodes[7] =
Point (.5,.5,-1.);
565 nodes[8] =
Point (0.,.5,-1.);
566 nodes[9] =
Point (0.,0.,0.);
567 nodes[10] =
Point (1.,0.,0.);
568 nodes[11] =
Point (0.,1.,0.);
569 nodes[12] =
Point (.5,0.,1.);
570 nodes[13] =
Point (.5,.5,1.);
571 nodes[14] =
Point (0.,.5,1.);
572 libmesh_fallthrough();
576 nodes[0] =
Point (0.,0.,-1.);
577 nodes[1] =
Point (1.,0.,-1.);
578 nodes[2] =
Point (0.,1.,-1.);
579 nodes[3] =
Point (0.,0.,1.);
580 nodes[4] =
Point (1.,0.,1.);
581 nodes[5] =
Point (0.,1.,1.);
592 libmesh_fallthrough();
597 nodes[13] =
Point (0.,0.,0.);
599 libmesh_fallthrough();
604 nodes[5] =
Point (0.,-1.,0.);
605 nodes[6] =
Point (1.,0.,0.);
606 nodes[7] =
Point (0.,1.,0.);
607 nodes[8] =
Point (-1,0.,0.);
610 nodes[9] =
Point (-.5,-.5,.5);
611 nodes[10] =
Point (.5,-.5,.5);
612 nodes[11] =
Point (.5,.5,.5);
613 nodes[12] =
Point (-.5,.5,.5);
615 libmesh_fallthrough();
620 nodes[0] =
Point (-1.,-1.,0.);
621 nodes[1] =
Point (1.,-1.,0.);
622 nodes[2] =
Point (1.,1.,0.);
623 nodes[3] =
Point (-1.,1.,0.);
625 nodes[4] =
Point (0.,0.,1.);
636 #ifdef LIBMESH_ENABLE_DEPRECATED 640 libmesh_deprecated();
642 libmesh_assert_greater_equal (eps, 0.);
644 const Real xi = p(0);
646 const Real eta = p(1);
651 const Real zeta = p(2);
653 const Real zeta = 0.;
660 return (!xi && !eta && !zeta);
667 if ((xi >= -1.-eps) &&
682 if ((xi >= 0.-eps) &&
684 ((xi + eta) <= 1.+eps))
699 if ((xi >= -1.-eps) &&
716 if ((xi >= 0.-eps) &&
719 ((xi + eta + zeta) <= 1.+eps))
741 if ((xi >= -1.-eps) &&
764 if ((xi >= 0.-eps) &&
768 ((xi + eta) <= 1.+eps))
788 if ((-eta - 1. + zeta <= 0.+eps) &&
789 ( xi - 1. + zeta <= 0.+eps) &&
790 ( eta - 1. + zeta <= 0.+eps) &&
791 ( -xi - 1. + zeta <= 0.+eps) &&
798 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 804 if ((xi >= -1.-eps) &&
820 if ((xi >= 0.-eps) &&
824 ((xi + eta) <= 1.+eps))
842 #endif // LIBMESH_ENABLE_DEPRECATED 861 os <<
"phi[i][j]: Shape function i at quadrature pt. j" << std::endl;
864 os <<
"dphi[i][j]: Shape function i's gradient at quadrature pt. j" << std::endl;
867 os <<
"XYZ locations of the quadrature pts." << std::endl;
870 os <<
"Values of JxW at the quadrature pts." << std::endl;
883 #ifdef LIBMESH_ENABLE_AMR 885 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 891 const unsigned int Dim = elem->
dim();
902 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 912 inf_fe_family_mapping_switch(2, inf_compute_node_constraints (constraints, elem) , ,;
break;);
917 inf_fe_family_mapping_switch(3, inf_compute_node_constraints (constraints, elem) , ,;
break;);
921 libmesh_error_msg(
"Invalid dim = " << Dim);
931 std::vector<const Node *> my_nodes, parent_nodes;
932 std::unique_ptr<const Elem> my_side, parent_side;
953 const unsigned int n_side_nodes = my_side->n_nodes();
956 my_nodes.reserve (n_side_nodes);
957 parent_nodes.clear();
958 parent_nodes.reserve (n_side_nodes);
960 for (
unsigned int n=0; n != n_side_nodes; ++n)
961 my_nodes.push_back(my_side->node_ptr(n));
963 for (
unsigned int n=0; n != n_side_nodes; ++n)
964 parent_nodes.push_back(parent_side->node_ptr(n));
966 for (
unsigned int my_side_n=0;
967 my_side_n < n_side_nodes;
977 const int side_max_order =
981 side_fe_type.
order = side_max_order;
989 const Node * my_node = my_nodes[my_side_n];
992 const Point & support_point = *my_node;
1000 for (
unsigned int their_side_n=0;
1001 their_side_n < n_side_nodes;
1009 parent_side.get()));
1011 const Node * their_node = parent_nodes[their_side_n];
1021 const Real their_mag = std::abs(their_value);
1025 if (their_mag > 0.999)
1027 libmesh_assert_equal_to (my_node, their_node);
1028 libmesh_assert_less (std::abs(their_value - 1.), 0.001);
1037 if (their_mag < 1.e-5)
1047 constraint_row.emplace(their_node, 0.);
1062 constraint_row.emplace(their_node, their_value);
1069 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS 1071 #endif // #ifdef LIBMESH_ENABLE_AMR 1075 #ifdef LIBMESH_ENABLE_PERIODIC 1077 #ifdef LIBMESH_ENABLE_NODE_CONSTRAINTS 1085 if (boundaries.empty())
1094 const unsigned int Dim = elem->
dim();
1100 std::vector<const Node *> my_nodes, neigh_nodes;
1101 std::unique_ptr<const Elem> my_side, neigh_side;
1105 std::vector<boundary_id_type> bc_ids;
1111 mesh.get_boundary_info().boundary_ids (elem, s, bc_ids);
1112 for (
const auto & boundary_id : bc_ids)
1120 unsigned int s_neigh;
1121 const Elem * neigh = boundaries.
neighbor(boundary_id, *point_locator, elem, s, &s_neigh);
1123 libmesh_error_msg_if
1124 (!neigh,
"PeriodicBoundaries can't find a periodic neighbor for element " <<
1125 elem->
id() <<
" side " << s);
1133 #ifdef LIBMESH_ENABLE_AMR 1135 #endif // #ifdef LIBMESH_ENABLE_AMR 1140 const unsigned int n_side_nodes = my_side->n_nodes();
1143 my_nodes.reserve (n_side_nodes);
1144 neigh_nodes.clear();
1145 neigh_nodes.reserve (n_side_nodes);
1147 for (
unsigned int n=0; n != n_side_nodes; ++n)
1148 my_nodes.push_back(my_side->node_ptr(n));
1150 for (
unsigned int n=0; n != n_side_nodes; ++n)
1151 neigh_nodes.push_back(neigh_side->node_ptr(n));
1157 std::vector<bool> skip_constraint(n_side_nodes,
false);
1159 for (
unsigned int my_side_n=0;
1160 my_side_n < n_side_nodes;
1166 const Node * my_node = my_nodes[my_side_n];
1174 if (constraints.count(my_node))
1176 skip_constraint[my_side_n] =
true;
1182 for (
unsigned int their_side_n=0;
1183 their_side_n < n_side_nodes;
1189 const Node * their_node = neigh_nodes[their_side_n];
1198 if (!constraints.count(their_node))
1202 constraints[their_node].first;
1204 for (
unsigned int orig_side_n=0;
1205 orig_side_n < n_side_nodes;
1211 const Node * orig_node = my_nodes[orig_side_n];
1213 if (their_constraint_row.count(orig_node))
1214 skip_constraint[orig_side_n] =
true;
1219 for (
unsigned int my_side_n=0;
1220 my_side_n < n_side_nodes;
1226 if (skip_constraint[my_side_n])
1229 const Node * my_node = my_nodes[my_side_n];
1235 const Point mapped_point =
1239 for (
unsigned int their_side_n=0;
1240 their_side_n < n_side_nodes;
1246 const Node * their_node = neigh_nodes[their_side_n];
1263 constraints[my_node].first;
1265 constraint_row.emplace(their_node, their_value);
1274 #endif // LIBMESH_ENABLE_NODE_CONSTRAINTS 1276 #endif // LIBMESH_ENABLE_PERIODIC
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
FEFamily family
The type of finite element.
ElemType
Defines an enum for geometric element types.
virtual void print_phi(std::ostream &os) const =0
Prints the value of each shape function at each quadrature point.
unsigned int _n_total_qp
The total number of quadrature points for the current configuration.
const Elem * parent() const
A Node is like a Point, but with more information.
static unsigned int n_dofs(const unsigned int dim, const FEType &fe_t, const ElemType t)
virtual Order default_side_order() const
const unsigned int invalid_uint
A number which is used quite often to represent an invalid or uninitialized value for an unsigned int...
bool shapes_on_quadrature
A flag indicating if current data structures correspond to quadrature rule points.
static void compute_node_constraints(NodeConstraints &constraints, const Elem *elem)
Computes the nodal constraint contributions (for non-conforming adapted meshes), using Lagrange geome...
IntRange< unsigned short > side_index_range() const
virtual std::unique_ptr< Elem > build_side_ptr(const unsigned int i)=0
static Point inverse_map(const unsigned int dim, const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true, const bool extra_checks=true)
We're using a class instead of a typedef to allow forward declarations and future flexibility...
This is the base class from which all geometric element types are derived.
PeriodicBoundaryBase * boundary(boundary_id_type id)
FEAbstract(const unsigned int dim, const FEType &fet)
Constructor.
std::ostream & operator<<(std::ostream &os, const OrderWrapper &order)
Overload stream operators.
The Node constraint storage format.
OrderWrapper order
The approximation order of the element.
The libMesh namespace provides an interface to certain functionality in the library.
static const unsigned int type_to_n_nodes_map[INVALID_ELEM]
This array maps the integer representation of the ElemType enum to the number of nodes in the element...
static unsigned int max_order(const FEType &fe_t, const ElemType &el_t)
static bool on_reference_element(const Point &p, const ElemType t, const Real eps=TOLERANCE)
This is the MeshBase class.
virtual Point get_corresponding_pos(const Point &pt) const =0
This function should be overridden by derived classes to define how one finds corresponding nodes on ...
virtual void print_dphi(std::ostream &os) const =0
Prints the value of each shape function's derivative at each quadrature point.
const dof_id_type n_nodes
const unsigned int dim
The dimensionality of the object.
void print_xyz(std::ostream &os) const
Prints the spatial location of each quadrature point (on the physical element).
static Real shape(const unsigned int dim, const FEType &fe_t, const ElemType t, const unsigned int i, const Point &p)
virtual unsigned int n_quadrature_points() const
static void get_refspace_nodes(const ElemType t, std::vector< Point > &nodes)
unsigned int n_points() const
QBase * qrule
A pointer to the quadrature rule employed.
This is the base class for point locators.
void print_JxW(std::ostream &os) const
Prints the Jacobian times the weight for each quadrature point.
std::string enum_to_string(const T e)
void print_info(std::ostream &os) const
Prints all the relevant information about the current element.
const Elem * neighbor_ptr(unsigned int i) const
unsigned int level() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual unsigned short dim() const =0
static std::unique_ptr< FEAbstract > build(const unsigned int dim, const FEType &type)
Builds a specific finite element type.
const Elem * neighbor(boundary_id_type boundary_id, const PointLocatorBase &point_locator, const Elem *e, unsigned int side, unsigned int *neigh_side=nullptr) const
virtual ~FEAbstract()
Destructor.
The base class for defining periodic boundaries.
std::map< const Node *, Real, std::less< const Node * >, Threads::scalable_allocator< std::pair< const Node *const, Real > > > NodeConstraintRow
A row of the Node constraint mapping.
virtual bool infinite() const =0
FEFamily
defines an enum for finite element families.
This class forms the foundation from which generic finite elements may be derived.
FEType fe_type
The finite element type for this object.
std::unique_ptr< FEMap > _fe_map
virtual Order default_order() const =0
Class contained in FE that encapsulates mapping (i.e.
A Point defines a location in LIBMESH_DIM dimensional Real space.
static void compute_periodic_node_constraints(NodeConstraints &constraints, const PeriodicBoundaries &boundaries, const MeshBase &mesh, const PointLocatorBase *point_locator, const Elem *elem)
Computes the node position constraint equation contributions (for meshes with periodic boundary condi...
static FEFamily map_fe_type(const Elem &elem)
spin_mutex spin_mtx
A convenient spin mutex object which can be used for obtaining locks.
const RemoteElem * remote_elem