LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Cubic Lagrangian finite elment on a triangular reference element. More...
#include <lf/uscalfe/uscalfe.h>
Public Member Functions | |
FeLagrangeO3Tria (const FeLagrangeO3Tria &)=default | |
FeLagrangeO3Tria (FeLagrangeO3Tria &&) noexcept=default | |
FeLagrangeO3Tria & | operator= (const FeLagrangeO3Tria &)=default |
FeLagrangeO3Tria & | operator= (FeLagrangeO3Tria &&) noexcept=default |
FeLagrangeO3Tria ()=default | |
~FeLagrangeO3Tria () override=default | |
base::RefEl | RefEl () const override |
Tells the type of reference cell underlying the parametric finite element. | |
unsigned | Degree () const override |
Quadratic Lagrangian finite elements rely on polynomials of degree 3. | |
size_type | NumRefShapeFunctions () const override |
Ten local shape functions are associated with a triangular cell in the case of cubic Lagrangian finite elements. | |
size_type | NumRefShapeFunctions (dim_t codim) const override |
One shape function attached to each node and two to each edge of the triangle. There is one interior shape function. | |
size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const override |
One shape function attached to each node and two to each edge of the triangle. There is one interior shape function. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | EvalReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override |
Point evaluation of reference shape functions. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | GradientsReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override |
Point evaluations of gradients of reference shape functions. | |
Eigen::MatrixXd | EvaluationNodes () const override |
Evaluation nodes are the three vertices of the triangle, two equispaced interio nodes for each edge, a the barycentre. | |
size_type | NumEvaluationNodes () const override |
Ten evaluation nodes. | |
![]() | |
virtual | ~ScalarReferenceFiniteElement ()=default |
dim_t | Dimension () const |
Returns the spatial dimension of the reference cell. | |
virtual Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > | NodalValuesToDofs (const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodvals) const |
Computes the linear combination of reference shape functions matching function values at evaluation nodes. | |
Additional Inherited Members | |
![]() | |
using | Scalar = SCALAR |
The underlying scalar type. | |
![]() | |
ScalarReferenceFiniteElement ()=default | |
ScalarReferenceFiniteElement (const ScalarReferenceFiniteElement &)=default | |
ScalarReferenceFiniteElement (ScalarReferenceFiniteElement &&) noexcept=default | |
ScalarReferenceFiniteElement & | operator= (const ScalarReferenceFiniteElement &)=default |
ScalarReferenceFiniteElement & | operator= (ScalarReferenceFiniteElement &&) noexcept=default |
![]() | |
template<class SCALAR > | |
void | PrintInfo (std::ostream &o, const ScalarReferenceFiniteElement< SCALAR > &srfe, unsigned int ctrl=0) |
Print information about a ScalarReferenceFiniteElement to the given stream. | |
template<typename SCALAR > | |
std::ostream & | operator<< (std::ostream &o, const ScalarReferenceFiniteElement< SCALAR > &fe_desc) |
Stream output operator: just calls the ScalarReferenceFiniteElement::print() method. | |
Cubic Lagrangian finite elment on a triangular reference element.
This is a specialization of ScalarReferenceFiniteElement. Refer to its documentation.
The 10 reference shape functions are combinations of the barycentric coordinate functions on the reference triangle.
The first three shape functions are associated with vertices, the next six with edges of the triangle. The last shape function is an interior shape function, see Example 2.6.1.7 for the location of interpolation nodes.
|
default |
|
defaultnoexcept |
|
default |
|
overridedefault |
|
inlineoverridevirtual |
Quadratic Lagrangian finite elements rely on polynomials of degree 3.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Point evaluation of reference shape functions.
The reference shape function are cardinal basis functions of the space of cubic polynomials with respect to the ten interpolation nodes supplied by the lattive drawn in Figure 111.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Evaluation nodes are the three vertices of the triangle, two equispaced interio nodes for each edge, a the barycentre.
The reference shape functions as implemented by this class are a cardinal basis of the space of cubic two-variate polynomials with respect to these evaluation nodes.
The numbering of evluation nodes is fixed by the numbering of the local shape functions.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Point evaluations of gradients of reference shape functions.
Compute values of gradients of reference shape functions at specified nodes
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Ten evaluation nodes.
Tell the number of evaluation (interpolation) nodes.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 1137 of file lagr_fe.h.
References lf::uscalfe::FeLagrangeO3Tria< SCALAR >::NumRefShapeFunctions().
|
inlineoverridevirtual |
Ten local shape functions are associated with a triangular cell in the case of cubic Lagrangian finite elements.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 969 of file lagr_fe.h.
Referenced by lf::uscalfe::FeLagrangeO3Tria< SCALAR >::NumEvaluationNodes().
|
inlineoverridevirtual |
One shape function attached to each node and two to each edge of the triangle. There is one interior shape function.
The number of interior reference shape functions for sub-entities of a particular co-dimension.
codim | co-dimension of the subentity |
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
One shape function attached to each node and two to each edge of the triangle. There is one interior shape function.
The number of interior reference shape functions for every sub-entity.
codim | do-dimension of the subentity |
subidx | local index of the sub-entity |
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
default |
|
defaultnoexcept |
|
inlineoverridevirtual |
Tells the type of reference cell underlying the parametric finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 955 of file lagr_fe.h.
References lf::base::RefEl::kTria().