24 #include "libmesh/fe_base.h" 25 #include "libmesh/int_range.h" 26 #include "libmesh/libmesh.h" 39 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 41 template <
unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
126 template <
unsigned int Dim, FEFamily T>
151 const unsigned int i,
166 const unsigned int i,
168 const bool add_p_level =
true);
182 const unsigned int i,
184 const bool add_p_level =
true);
198 const unsigned int i,
199 const std::vector<Point> & p,
200 std::vector<OutputShape> & v,
201 const bool add_p_level =
true);
216 const std::vector<Point> & p,
217 std::vector<std::vector<OutputShape> > & v,
218 const bool add_p_level =
true);
229 const unsigned int i,
230 const unsigned int j,
243 const unsigned int i,
244 const unsigned int j,
246 const bool add_p_level =
true);
259 const unsigned int i,
260 const unsigned int j,
262 const bool add_p_level =
true);
276 const unsigned int i,
277 const unsigned int j,
278 const std::vector<Point> & p,
279 std::vector<OutputShape> & v,
280 const bool add_p_level =
true);
295 const std::vector<Point> & p,
296 std::vector<std::vector<OutputShape>> * comps[3],
297 const bool add_p_level =
true);
299 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 322 const unsigned int i,
323 const unsigned int j,
350 const unsigned int i,
351 const unsigned int j,
353 const bool add_p_level =
true);
377 const unsigned int i,
378 const unsigned int j,
380 const bool add_p_level =
true);
382 #endif //LIBMESH_ENABLE_SECOND_DERIVATIVES 390 const std::vector<Number> & elem_soln,
392 bool add_p_level =
true,
393 const unsigned vdim = 1);
402 const unsigned int side,
403 const std::vector<Number> & elem_soln,
404 std::vector<Number> & nodal_soln_on_side,
405 bool add_p_level =
true,
406 const unsigned vdim = 1);
467 const unsigned int n);
477 const unsigned int n);
520 std::vector<unsigned int> & di,
521 bool add_p_level=
true);
531 std::vector<unsigned int> & di,
532 bool add_p_level=
true);
537 const bool secure =
true)
544 const std::vector<Point> & physical_points,
545 std::vector<Point> & reference_points,
547 const bool secure =
true)
551 tolerance, secure, secure);
565 const std::vector<Point> *
const pts =
nullptr,
566 const std::vector<Real> *
const weights =
nullptr)
override;
573 const std::vector<Point> & pts,
574 const std::vector<Real> & JxW)
override;
592 const unsigned int side,
594 const std::vector<Point> *
const pts =
nullptr,
595 const std::vector<Real> *
const weights =
nullptr)
override;
607 const unsigned int edge,
609 const std::vector<Point> *
const pts =
nullptr,
610 const std::vector<Real> *
const weights =
nullptr)
override;
618 const unsigned int s,
619 const std::vector<Point> & reference_side_points,
620 std::vector<Point> & reference_points)
override;
628 const unsigned int e,
629 const std::vector<Point> & reference_edge_points,
630 std::vector<Point> & reference_points);
639 #ifdef LIBMESH_ENABLE_AMR 648 const unsigned int variable_number,
650 #endif // #ifdef LIBMESH_ENABLE_AMR 661 const Point & reference_point)
664 return FEMap::map(Dim, elem, reference_point);
668 const Point & reference_point)
675 const Point & reference_point)
682 const Point & reference_point)
688 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 693 template <
unsigned int friend_Dim, FEFamily friend_T_radial, InfMapType friend_T_map>
714 const std::vector<Point> & p,
715 std::vector<std::vector<OutputShape>> * comps[3],
716 const bool add_p_level =
true);
724 #ifdef LIBMESH_ENABLE_INFINITE_ELEMENTS 731 const Elem * e)
override;
740 const unsigned int i,
741 const std::vector<Point> & p,
742 std::vector<OutputShape> & v,
743 const bool add_p_level =
true)
745 libmesh_assert_equal_to(p.size(), v.size());
755 const std::vector<Point> & p,
756 std::vector<std::vector<OutputShape>> & v,
757 const bool add_p_level =
true)
761 libmesh_assert_equal_to ( p.size(), v[i].size() );
771 const unsigned int i,
772 const unsigned int j,
773 const std::vector<Point> & p,
774 std::vector<OutputShape> & v,
775 const bool add_p_level =
true)
777 libmesh_assert_equal_to(p.size(), v.size());
786 const unsigned int side,
787 const std::vector<Number> & elem_soln,
788 std::vector<Number> & nodal_soln_on_side,
789 bool add_p_level =
true,
790 const unsigned vdim = 1);
828 template <
unsigned int Dim>
852 template <
unsigned int Dim>
866 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 890 const unsigned int i,
892 const bool add_p_level);
897 const unsigned int i,
898 const unsigned int j,
900 const bool add_p_level);
902 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 906 const unsigned int i,
907 const unsigned int j,
909 const bool add_p_level);
935 const std::vector<Point> *
const pts =
nullptr,
936 const std::vector<Real> *
const weights =
nullptr)
override;
945 const std::vector<Point> *
const =
nullptr,
946 const std::vector<Real> *
const =
nullptr)
override 947 { libmesh_not_implemented(); }
964 const Elem * elem)
override;
983 const unsigned int j,
987 #ifdef LIBMESH_ENABLE_SECOND_DERIVATIVES 995 const unsigned int j,
1000 #endif // LIBMESH_ENABLE_SECOND_DERIVATIVE 1011 const unsigned int valence);
1019 unsigned int valence);
1031 template <
unsigned int Dim>
1055 template <
unsigned int Dim>
1079 template <
unsigned int Dim>
1098 template <
unsigned int Dim>
1121 template <
unsigned int Dim>
1140 template <
unsigned int Dim>
1166 template <
unsigned int Dim>
1186 const std::vector<Point> *
const pts =
nullptr,
1187 const std::vector<Real> *
const weights =
nullptr)
override 1195 const unsigned int side,
1197 const std::vector<Point> *
const pts =
nullptr,
1198 const std::vector<Real> *
const weights =
nullptr)
override;
1211 const Elem * e)
override;
1230 const std::vector<Real> & weights);
1242 template <
unsigned int Dim>
1265 template <
unsigned int Dim>
1289 template <
unsigned int Dim>
1312 template <
unsigned int Dim>
1335 template <
unsigned int Dim>
1356 template <
unsigned int Dim>
1378 template <
unsigned int Dim>
1401 template <
unsigned int Dim>
1419 namespace FiniteElements
1527 template <
typename OutputShape>
1530 const unsigned int i,
1531 const unsigned int j,
1533 const bool add_p_level,
1534 OutputShape(*shape_func)
1536 const unsigned int,
const Point &,
1539 template <
typename OutputShape>
1542 const unsigned int i,
1543 const unsigned int j,
1545 OutputShape(*shape_func)
1547 const unsigned int,
const Point &));
1549 template <
typename OutputShape>
1553 const unsigned int i,
1554 const unsigned int j,
1556 OutputShape(*shape_func)
1558 const Elem *,
const unsigned int,
1562 template <
typename OutputShape>
1566 const unsigned int i,
1567 const unsigned int j,
1569 const bool add_p_level,
1570 OutputShape(*deriv_func)
1572 const unsigned int,
const unsigned int,
1573 const Point &,
const bool));
1575 template <
typename OutputShape>
1578 const unsigned int i,
1579 const unsigned int j,
1581 OutputShape(*deriv_func)
1587 template <
typename OutputShape>
1591 const unsigned int i,
1592 const unsigned int j,
1594 OutputShape(*deriv_func)
1606 const std::vector<Number> & elem_soln,
1607 std::vector<Number> & nodal_soln,
1608 bool add_p_level =
true);
1622 const FEType underlying_fe_type,
1623 std::vector<std::vector<Real>> & shapes,
1624 const std::vector<Point> & p,
1625 const bool add_p_level);
1631 std::vector<std::vector<Real>> & shapes,
1632 std::vector<std::vector<std::vector<Real>>> & derivs,
1633 const std::vector<Point> & p,
1634 const bool add_p_level);
1637 const FEType underlying_fe_type,
1638 const unsigned int i,
1640 const bool add_p_level);
1643 const FEType underlying_fe_type,
1644 const unsigned int i,
1645 const unsigned int j,
1647 const bool add_p_level);
1650 const FEType underlying_fe_type,
1651 const unsigned int i,
1652 const unsigned int j,
1654 const bool add_p_level);
1657 const FEType underlying_fe_type,
1658 const std::vector<Point> & p,
1659 std::vector<std::vector<Real>> & v,
1660 const bool add_p_level);
1662 template <
typename OutputShape>
1664 const FEType underlying_fe_type,
1665 const std::vector<Point> & p,
1666 std::vector<std::vector<OutputShape>> * comps[3],
1667 const bool add_p_level);
1674 #define LIBMESH_DEFAULT_NDOFS(MyType) \ 1675 template <> unsigned int FE<0,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \ 1676 template <> unsigned int FE<1,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \ 1677 template <> unsigned int FE<2,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \ 1678 template <> unsigned int FE<3,MyType>::n_dofs(const ElemType t, const Order o) { return MyType##_n_dofs(t, o); } \ 1680 template <> unsigned int FE<0,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \ 1681 template <> unsigned int FE<1,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \ 1682 template <> unsigned int FE<2,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \ 1683 template <> unsigned int FE<3,MyType>::n_dofs(const Elem * e, const Order o) { return MyType##_n_dofs(e, o); } \ 1685 template <> unsigned int FE<0,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \ 1686 template <> unsigned int FE<1,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \ 1687 template <> unsigned int FE<2,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \ 1688 template <> unsigned int FE<3,MyType>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(t, o, n); } \ 1690 template <> unsigned int FE<0,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \ 1691 template <> unsigned int FE<1,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \ 1692 template <> unsigned int FE<2,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \ 1693 template <> unsigned int FE<3,MyType>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return MyType##_n_dofs_at_node(e, o, n); } \ 1695 template <> unsigned int FE<0,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \ 1696 template <> unsigned int FE<1,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \ 1697 template <> unsigned int FE<2,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \ 1698 template <> unsigned int FE<3,MyType>::n_dofs_per_elem(const ElemType t, const Order o) { return MyType##_n_dofs_per_elem(t, o); } \ 1700 template <> unsigned int FE<0,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \ 1701 template <> unsigned int FE<1,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \ 1702 template <> unsigned int FE<2,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } \ 1703 template <> unsigned int FE<3,MyType>::n_dofs_per_elem(const Elem & e, const Order o) { return MyType##_n_dofs_per_elem(e, o); } 1706 #define LIBMESH_DEFAULT_VEC_NDOFS(MyType) \ 1707 template <> unsigned int FE<0,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return FE<0,MyType>::n_dofs(t, o); } \ 1708 template <> unsigned int FE<1,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return FE<1,MyType>::n_dofs(t, o); } \ 1709 template <> unsigned int FE<2,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return 2*FE<2,MyType>::n_dofs(t, o); } \ 1710 template <> unsigned int FE<3,MyType##_VEC>::n_dofs(const ElemType t, const Order o) { return 3*FE<3,MyType>::n_dofs(t, o); } \ 1712 template <> unsigned int FE<0,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return FE<0,MyType>::n_dofs(e, o); } \ 1713 template <> unsigned int FE<1,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return FE<1,MyType>::n_dofs(e, o); } \ 1714 template <> unsigned int FE<2,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return 2*FE<2,MyType>::n_dofs(e, o); } \ 1715 template <> unsigned int FE<3,MyType##_VEC>::n_dofs(const Elem * e, const Order o) { return 3*FE<3,MyType>::n_dofs(e, o); } \ 1717 template <> unsigned int FE<0,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<0,MyType>::n_dofs_at_node(t, o, n); } \ 1718 template <> unsigned int FE<1,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return FE<1,MyType>::n_dofs_at_node(t, o, n); } \ 1719 template <> unsigned int FE<2,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return 2*FE<2,MyType>::n_dofs_at_node(t, o, n); } \ 1720 template <> unsigned int FE<3,MyType##_VEC>::n_dofs_at_node(const ElemType t, const Order o, const unsigned int n) { return 3*FE<3,MyType>::n_dofs_at_node(t, o, n); } \ 1722 template <> unsigned int FE<0,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<0,MyType>::n_dofs_at_node(e.type(), o, n); } \ 1723 template <> unsigned int FE<1,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return FE<1,MyType>::n_dofs_at_node(e.type(), o, n); } \ 1724 template <> unsigned int FE<2,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return 2*FE<2,MyType>::n_dofs_at_node(e.type(), o, n); } \ 1725 template <> unsigned int FE<3,MyType##_VEC>::n_dofs_at_node(const Elem & e, const Order o, const unsigned int n) { return 3*FE<3,MyType>::n_dofs_at_node(e.type(), o, n); } \ 1727 template <> unsigned int FE<0,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<0,MyType>::n_dofs_per_elem(t, o); } \ 1728 template <> unsigned int FE<1,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return FE<1,MyType>::n_dofs_per_elem(t, o); } \ 1729 template <> unsigned int FE<2,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 2*FE<2,MyType>::n_dofs_per_elem(t, o); } \ 1730 template <> unsigned int FE<3,MyType##_VEC>::n_dofs_per_elem(const ElemType t, const Order o) { return 3*FE<3,MyType>::n_dofs_per_elem(t, o); } \ 1732 template <> unsigned int FE<0,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<0,MyType>::n_dofs_per_elem(e.type(), o); } \ 1733 template <> unsigned int FE<1,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return FE<1,MyType>::n_dofs_per_elem(e.type(), o); } \ 1734 template <> unsigned int FE<2,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 2*FE<2,MyType>::n_dofs_per_elem(e.type(), o); } \ 1735 template <> unsigned int FE<3,MyType##_VEC>::n_dofs_per_elem(const Elem & e, const Order o) { return 3*FE<3,MyType>::n_dofs_per_elem(e.type(), o); } 1738 #define LIBMESH_DEFAULT_VECTORIZED_FE(MyDim, MyType) \ 1740 void FE<MyDim,MyType>::all_shapes \ 1741 (const Elem * elem, \ 1743 const std::vector<Point> & p, \ 1744 std::vector<std::vector<OutputShape>> & v, \ 1745 const bool add_p_level) \ 1747 FE<MyDim,MyType>::default_all_shapes \ 1748 (elem,o,p,v,add_p_level); \ 1752 void FE<MyDim,MyType>::shapes \ 1753 (const Elem * elem, \ 1755 const unsigned int i, \ 1756 const std::vector<Point> & p, \ 1757 std::vector<OutputShape> & v, \ 1758 const bool add_p_level) \ 1760 FE<MyDim,MyType>::default_shapes \ 1761 (elem,o,i,p,v,add_p_level); \ 1765 void FE<MyDim,MyType>::shape_derivs \ 1766 (const Elem * elem, \ 1768 const unsigned int i, \ 1769 const unsigned int j, \ 1770 const std::vector<Point> & p, \ 1771 std::vector<OutputShape> & v, \ 1772 const bool add_p_level) \ 1774 FE<MyDim,MyType>::default_shape_derivs \ 1775 (elem,o,i,j,p,v,add_p_level); \ 1779 void FE<MyDim,MyType>::all_shape_derivs \ 1780 (const Elem * elem, \ 1782 const std::vector<Point> & p, \ 1783 std::vector<std::vector<OutputShape>> * comps[3], \ 1784 const bool add_p_level) \ 1786 FE<MyDim,MyType>::default_all_shape_derivs \ 1787 (elem,o,p,comps,add_p_level); \ 1791 #endif // LIBMESH_FE_H
FELagrangeVec objects are used for working with vector-valued finite elements.
class FEType hides (possibly multiple) FEFamily and approximation orders, thereby enabling specialize...
FELagrange(const FEType &fet)
Constructor.
static Point map_xi(const Elem *elem, const Point &reference_point)
FESubdivision(const FEType &fet)
Constructor.
FEMonomial(const FEType &fet)
Constructor.
static unsigned int n_dofs(const ElemType t, const Order o)
ElemType
Defines an enum for geometric element types.
static Real regular_shape(const unsigned int i, const Real v, const Real w)
Order
defines an enum for polynomial orders.
static void dofs_on_edge(const Elem *const elem, const Order o, unsigned int e, std::vector< unsigned int > &di, bool add_p_level=true)
Fills the vector di with the local degree of freedom indices associated with edge e of element elem...
A specific instantiation of the FEBase class.
FENedelecOne objects are used for working with vector-valued Nedelec finite elements of the first kin...
virtual void side_map(const Elem *elem, const Elem *side, const unsigned int s, const std::vector< Point > &reference_side_points, std::vector< Point > &reference_points) override
Computes the reference space quadrature points on the side of an element based on the side quadrature...
static OutputShape shape(const ElemType t, const Order o, const unsigned int i, const Point &p)
FE< 1, MONOMIAL > FEMonomial1D
Convenient definition for a 1D Monomial finite element.
FE< 2, L2_HIERARCHIC > FEL2Hierarchic2D
Convenient definition for a 2D Discontinuous Hierarchic finite element.
virtual void attach_quadrature_rule(QBase *q) override
Provides the class with the quadrature rule, which provides the locations (on a reference element) wh...
FE< 3, L2_HIERARCHIC > FEL2Hierarchic3D
Convenient definition for a 3D Discontinuous Hierarchic finite element.
static void init_subdivision_matrix(DenseMatrix< Real > &A, unsigned int valence)
Builds the subdivision matrix A for the Loop scheme.
void compute_face_values(const Elem *elem, const Elem *side, const std::vector< Real > &weights)
Compute the map & shape functions for this face.
void rational_fe_weighted_shapes_derivs(const Elem *elem, const FEType fe_type, std::vector< std::vector< Real >> &shapes, std::vector< std::vector< std::vector< Real >>> &derivs, const std::vector< Point > &p, const bool add_p_level)
Lagrange finite elements.
Monomial finite elements.
void init_dual_shape_functions(unsigned int n_shapes, unsigned int n_qp)
Init dual_phi and potentially dual_dphi, dual_d2phi.
FEClough(const FEType &fet)
Constructor.
static void loop_subdivision_mask(std::vector< Real > &weights, const unsigned int valence)
Fills the vector weights with the weight coefficients of the Loop subdivision mask for evaluating the...
static constexpr Real TOLERANCE
std::vector< bool > cached_faces
FE< 1, L2_LAGRANGE > FEL2Lagrange1D
Convenient definition for a 1D Discontinuous Lagrange finite element.
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)
static void all_shape_derivs(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> *comps[3], const bool add_p_level=true)
Fills comps with dphidxi (and in higher dimensions, eta/zeta) derivative component values for all sha...
static void default_all_shape_derivs(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> *comps[3], const bool add_p_level=true)
A default implementation for all_shape_derivs.
FEHierarchicVec objects are used for working with vector-valued high-order piecewise-continuous finit...
FERaviartThomas objects are used for working with vector-valued Raviart-Thomas finite elements...
static Point map(const Elem *elem, const Point &reference_point)
FERaviartThomas(const FEType &fet)
Constructor.
This is the base class from which all geometric element types are derived.
FEGenericBase< typename FEOutputType< T >::type >::OutputShape OutputShape
FE(const FEType &fet)
Constructor.
static unsigned int n_dofs_at_node(const ElemType t, const Order o, const unsigned int n)
static OutputShape shape_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
FE< 3, L2_LAGRANGE > FEL2Lagrange3D
Convenient definition for a 3D Discontinuous Lagrange finite element.
FEL2RaviartThomas(const FEType &fet)
Constructor.
virtual bool shapes_need_reinit() const override
FE< 1, HIERARCHIC > FEHierarchic1D
Convenient definition for a 1D Hierarchic finite element.
The libMesh namespace provides an interface to certain functionality in the library.
static Real hermite_raw_shape_second_deriv(const unsigned int basis_num, const Real xi)
1D hermite functions on unit interval
FEL2RaviartThomas objects are used for working with vector-valued discontinuous Raviart-Thomas finite...
void cache(const Elem *elem)
Repopulate the element cache with the node locations, edge and face orientations of the element elem...
FEXYZ(const FEType &fet)
Constructor.
FE< 2, LAGRANGE > FELagrange2D
Convenient definition for a 2D Lagrange finite element.
FELagrangeVec(const FEType &fet)
Constructor.
static void default_all_shapes(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> &v, const bool add_p_level=true)
A default implementation for all_shapes.
virtual void compute_shape_functions(const Elem *elem, const std::vector< Point > &qp) override
After having updated the jacobian and the transformation from local to global coordinates in FEAbstra...
FEL2LagrangeVec(const FEType &fet)
Constructor.
static Real hermite_raw_shape(const unsigned int basis_num, const Real xi)
virtual void attach_quadrature_rule(QBase *q) override
Provides the class with the quadrature rule, which provides the locations (on a reference element) wh...
FE< 3, LAGRANGE > FELagrange3D
Convenient definition for a 3D Lagrange finite element.
FEClough< 2 > FEClough2D
Convenient definition for a 2D Clough-Tocher finite element.
FEL2Lagrange(const FEType &fet)
Constructor.
FEHierarchicVec objects are used for working with vector-valued high-order finite elements...
virtual bool is_hierarchic() const override
OutputShape fe_fdm_deriv(const Elem *elem, const Order order, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level, OutputShape(*shape_func)(const Elem *, const Order, const unsigned int, const Point &, const bool))
Helper functions for finite differenced derivatives in cases where analytical calculations haven't be...
This class handles the numbering of degrees of freedom on a mesh.
FEMonomialVec(const FEType &fet)
Constructor.
FEHierarchic(const FEType &fet)
Constructor.
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
This is at the core of this class.
static Point map_eta(const Elem *elem, const Point &reference_point)
virtual unsigned int n_shape_functions() const override
A specific instantiation of the FEBase class.
static void dofs_on_side(const Elem *const elem, const Order o, unsigned int s, std::vector< unsigned int > &di, bool add_p_level=true)
Fills the vector di with the local degree of freedom indices associated with side s of element elem...
static Real regular_shape_second_deriv(const unsigned int i, const unsigned int j, const Real v, const Real w)
FEHermite(const FEType &fet)
Constructor.
FE< 3, HIERARCHIC > FEHierarchic3D
Convenient definition for a 3D Hierarchic finite element.
Discontinuous Hierarchic finite elements.
void rational_all_shape_derivs(const Elem &elem, const FEType underlying_fe_type, const std::vector< Point > &p, std::vector< std::vector< OutputShape >> *comps[3], const bool add_p_level)
unsigned int monomial_n_dofs(const ElemType t, const Order o)
Helper functions for Discontinuous-Pn type basis functions.
FEL2Hierarchic(const FEType &fet)
Constructor.
FE< 2, MONOMIAL > FEMonomial2D
Convenient definition for a 2D Monomial finite element.
static Real regular_shape_deriv(const unsigned int i, const unsigned int j, const Real v, const Real w)
FEOutputType< T >::type OutputShape
Convenient typedefs for gradients of output, hessians of output, and potentially-complex-valued versi...
OutputShape fe_fdm_second_deriv(const Elem *elem, const Order order, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level, OutputShape(*deriv_func)(const Elem *, const Order, const unsigned int, const unsigned int, const Point &, const bool))
static void default_shape_derivs(const Elem *elem, const Order o, const unsigned int i, const unsigned int j, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
A default implementation for shape_derivs.
static unsigned int n_dofs_per_elem(const ElemType t, const Order o)
virtual FEContinuity get_continuity() const override
FENedelecOne(const FEType &fet)
Constructor.
void lagrange_nodal_soln(const Elem *elem, const Order order, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, bool add_p_level=true)
Helper functions for Lagrange-based basis functions.
FEHierarchicVec(const FEType &fet)
Constructor.
std::vector< bool > cached_edges
static void side_nodal_soln(const Elem *elem, const Order o, const unsigned int side, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln_on_side, bool add_p_level=true, const unsigned vdim=1)
Build the nodal soln on one side from the (full) element soln.
Discontinuous Lagrange finite elements.
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...
static void nodal_soln(const Elem *elem, const Order o, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln, bool add_p_level=true, const unsigned vdim=1)
Build the nodal soln from the element soln.
Real rational_fe_shape(const Elem &elem, const FEType underlying_fe_type, const unsigned int i, const Point &p, const bool add_p_level)
Real rational_fe_shape_second_deriv(const Elem &elem, const FEType underlying_fe_type, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
static unsigned int n_shape_functions(const ElemType t, const Order o)
bool matches_cache(const Elem *elem)
Check if the node locations, edge and face orientations held in the element cache match those of elem...
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual void edge_map(const Elem *elem, const Elem *edge, const unsigned int e, const std::vector< Point > &reference_edge_points, std::vector< Point > &reference_points)
Computes the reference space quadrature points on the side of an element based on the edge quadrature...
FE< 2, HIERARCHIC > FEHierarchic2D
Convenient definition for a 2D Hierarchic finite element.
static Point map_zeta(const Elem *elem, const Point &reference_point)
Hierarchic finite elements.
FEContinuity
defines an enum for finite element types to libmesh_assert a certain level (or type? Hcurl?) of continuity.
static Point map(const unsigned int dim, const Elem *elem, const Point &reference_point)
static void inverse_map(const Elem *elem, const std::vector< Point > &physical_points, std::vector< Point > &reference_points, const Real tolerance=TOLERANCE, const bool secure=true)
FE< 3, MONOMIAL > FEMonomial3D
Convenient definition for a 3D Monomial finite element.
virtual void init_shape_functions(const std::vector< Point > &qp, const Elem *e)
Update the various member data fields phi, dphidxi, dphideta, dphidzeta, etc.
static void shape_derivs(const Elem *elem, const Order o, const unsigned int i, const unsigned int j, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
Fills v with the derivative of the shape function, evaluated at all points p.
FEMonomialVec objects are used for working with vector-valued discontinuous finite elements...
static Point map_deriv(const unsigned int dim, const Elem *elem, const unsigned int j, const Point &reference_point)
static Real hermite_raw_shape_deriv(const unsigned int basis_num, const Real xi)
FEL2LagrangeVec objects are used for working with vector-valued finite elements.
FE< 2, L2_LAGRANGE > FEL2Lagrange2D
Convenient definition for a 2D Discontinuous Lagrange finite element.
FE< 1, LAGRANGE > FELagrange1D
Convenient definition for a 1D Lagrange finite element.
static void shapes(const Elem *elem, const Order o, const unsigned int i, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
Fills v with the values of the shape function, evaluated at all points p.
The FEScalar class is used for working with SCALAR variables.
virtual void reinit_default_dual_shape_coeffs(const Elem *elem) override
This computes the default dual shape function coefficients.
FEScalar(const FEType &fet)
Constructor.
virtual void reinit_dual_shape_coeffs(const Elem *elem, const std::vector< Point > &pts, const std::vector< Real > &JxW) override
This re-computes the dual shape function coefficients.
FE< 1, L2_HIERARCHIC > FEL2Hierarchic1D
Convenient definition for a 1D Discontinuous Hierarchic finite element.
FEL2HierarchicVec(const FEType &fet)
Constructor.
void rational_all_shapes(const Elem &elem, const FEType underlying_fe_type, const std::vector< Point > &p, std::vector< std::vector< Real >> &v, const bool add_p_level)
std::vector< Point > cached_nodes
Vectors holding the node locations, edge and face orientations of the last element we cached...
virtual void init_shape_functions(const std::vector< Point > &qp, const Elem *elem) override
Update the various member data fields phi, dphidxi, dphideta, dphidzeta, etc.
Clough-Tocher finite elements.
ElemType last_side
The last side and last edge we did a reinit on.
static Point inverse_map(const Elem *elem, const Point &p, const Real tolerance=TOLERANCE, const bool secure=true)
void rational_fe_weighted_shapes(const Elem *elem, const FEType underlying_fe_type, std::vector< std::vector< Real >> &shapes, const std::vector< Point > &p, const bool add_p_level)
Helper functions for rational basis functions.
virtual void edge_reinit(const Elem *elem, const unsigned int edge, const Real tolerance=TOLERANCE, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
Reinitializes all the physical element-dependent data based on the edge.
virtual void init_shape_functions(const std::vector< Point > &qp, const Elem *e) override
Update the various member data fields phi, dphidxi, dphideta, dphidzeta, etc.
A Point defines a location in LIBMESH_DIM dimensional Real space.
The constraint matrix storage format.
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
This is at the core of this class.
virtual void reinit(const Elem *elem, const std::vector< Point > *const pts=nullptr, const std::vector< Real > *const weights=nullptr) override
Explicitly call base class method.
The QBase class provides the basic functionality from which various quadrature rules can be derived...
auto index_range(const T &sizable)
Helper function that returns an IntRange<std::size_t> representing all the indices of the passed-in v...
This class forms the foundation from which generic finite elements may be derived.
static OutputShape shape_second_deriv(const ElemType t, const Order o, const unsigned int i, const unsigned int j, const Point &p)
Most finite element types in libMesh are scalar-valued.
static void default_shapes(const Elem *elem, const Order o, const unsigned int i, const std::vector< Point > &p, std::vector< OutputShape > &v, const bool add_p_level=true)
A default implementation for shapes.
static void default_side_nodal_soln(const Elem *elem, const Order o, const unsigned int side, const std::vector< Number > &elem_soln, std::vector< Number > &nodal_soln_on_side, bool add_p_level=true, const unsigned vdim=1)
A default implementation for side_nodal_soln.
Real rational_fe_shape_deriv(const Elem &elem, const FEType underlying_fe_type, const unsigned int i, const unsigned int j, const Point &p, const bool add_p_level)
virtual void init_base_shape_functions(const std::vector< Point > &qp, const Elem *e) override
Initialize the data fields for the base of an an infinite element.
static void all_shapes(const Elem *elem, const Order o, const std::vector< Point > &p, std::vector< std::vector< OutputShape > > &v, const bool add_p_level=true)
Fills v[i][qp] with the values of the shape functions, evaluated at all points in p...
virtual void reinit(const Elem *, const unsigned int, const Real=TOLERANCE, const std::vector< Point > *const =nullptr, const std::vector< Real > *const =nullptr) override
This prevents some compilers being confused by partially overriding this virtual function.