LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Interface class for parametric scalar valued finite elements. More...
#include <lf/fe/fe.h>
Public Types | |
using | Scalar = SCALAR |
The underlying scalar type. | |
Public Member Functions | |
virtual | ~ScalarReferenceFiniteElement ()=default |
virtual base::RefEl | RefEl () const =0 |
Tells the type of reference cell underlying the parametric finite element. | |
virtual unsigned int | Degree () const =0 |
Request the maximal polynomial degree of the basis functions in this finite element. | |
dim_t | Dimension () const |
Returns the spatial dimension of the reference cell. | |
virtual size_type | NumRefShapeFunctions () const |
Total number of reference shape functions associated with the reference cell. | |
virtual size_type | NumRefShapeFunctions (dim_t codim) const |
The number of interior reference shape functions for sub-entities of a particular co-dimension. | |
virtual size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const =0 |
The number of interior reference shape functions for every sub-entity. | |
virtual Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | EvalReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const =0 |
Evaluation of all reference shape functions in a number of points. | |
virtual Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | GradientsReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const =0 |
Computation of the gradients of all reference shape functions in a number of points. | |
virtual Eigen::MatrixXd | EvaluationNodes () const =0 |
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom, nodal interpolation in the simplest case. | |
virtual size_type | NumEvaluationNodes () const =0 |
Tell the number of evaluation (interpolation) nodes. | |
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. | |
Protected Member Functions | |
ScalarReferenceFiniteElement ()=default | |
ScalarReferenceFiniteElement (const ScalarReferenceFiniteElement &)=default | |
ScalarReferenceFiniteElement (ScalarReferenceFiniteElement &&) noexcept=default | |
ScalarReferenceFiniteElement & | operator= (const ScalarReferenceFiniteElement &)=default |
ScalarReferenceFiniteElement & | operator= (ScalarReferenceFiniteElement &&) noexcept=default |
Related Symbols | |
(Note that these are not member symbols.) | |
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. | |
Interface class for parametric scalar valued finite elements.
SCALAR | underlying scalar type, usually either double or complex<double> |
A scalar parametric finite element is defined through a set of reference shape functions (RSFs) on a particular reference entity, cf. definition 2.8.3.1.
Each reference shape function is associated with a unique sub-entity of the reference entity according to Equation 2.8.1.6.
Specializations of this class support the evaluation of RSFs in arbitrary points in the reference element and the computation of their gradients. The also provide local components for the defintion of nodal interpolants.
This class is discussed in detail in Paragraph 2.8.3.25.
The numbering of reference shape functions is done according to the following convention, see also Paragraph 2.7.4.11
The numbering scheme must the same as that adopted for the definition of local-to-global maps, see the lf::assemble::DofHandler class.
For triangular cubic Lagrangian finite elements there is a single reference shape function associated with each vertex, two reference shape functions belonging to every edge and a single interior reference shape function. In this case the RSFs are numbered a follows
The rules governing this numbering in LehrFEM++ are explained above and in Paragraph 2.7.4.11.
Definition at line 82 of file scalar_reference_finite_element.h.
using lf::fe::ScalarReferenceFiniteElement< SCALAR >::Scalar = SCALAR |
The underlying scalar type.
Definition at line 98 of file scalar_reference_finite_element.h.
|
protecteddefault |
|
protecteddefault |
|
protecteddefaultnoexcept |
|
virtualdefault |
Virtual destructor
|
pure virtual |
Request the maximal polynomial degree of the basis functions in this finite element.
Implemented in lf::fe::FePoint< SCALAR >, lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Quad< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, lf::uscalfe::FeLagrangeO3Quad< SCALAR >, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >.
Referenced by lf::fe::ScalarReferenceFiniteElement< SCALAR >::PrintInfo().
|
inline |
Returns the spatial dimension of the reference cell.
Definition at line 119 of file scalar_reference_finite_element.h.
References lf::base::RefEl::Dimension(), and lf::fe::ScalarReferenceFiniteElement< SCALAR >::RefEl().
Referenced by lf::fe::ScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions().
|
pure virtual |
Evaluation of all reference shape functions in a number of points.
refcoords | coordinates of N points in the reference cell passed as columns of a matrix of size dim x N, where dim is the dimension of the reference element, that is =0 for points, =1 for edges, =2 for triangles and quadrilaterals |
NumRefShapeFunctions() x refcoords.cols()
which contains the shape functions evaluated at every quadrature point.Concerning the numbering of local shape functions, please consult the documentation of lf::assemble::DofHandler or the documentation of the class.
There are three reference shape functions \(\hat{b}^1,\hat{b}^2,\hat{b}^3\) associated with the vertices of the reference triangle. Let us assume that the refcoords
argument is a 2x2 matrix \([\mathbf{x}_1\;\mathbf{x}_2]\), which corresponds to passing the coordinates of two points in the reference triangle. Then this method will return a 3x2
matrix:
\[ \begin{pmatrix}\hat{b}^1(\mathbf{x}_1) & \hat{b}^1(\mathbf{x}_2) \\ \hat{b}^2(\mathbf{x}_1) & \hat{b}^2(\mathbf{x}_2) \\ \hat{b}^3(\mathbf{x}_1)\ & \hat{b}^3(\mathbf{x}_2) \end{pmatrix} \]
Implemented in lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >, lf::fe::FePoint< SCALAR >, lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Quad< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, and lf::uscalfe::FeLagrangeO3Quad< SCALAR >.
|
pure virtual |
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom, nodal interpolation in the simplest case.
Every parametric scalar finite element implicitly defines a local interpolation operator by duality with the reference shape functions. This interpolation operator can be realized through evaluations at certain evaluation nodes, which are provided by this method.
The evaluation points must satisfy the following requirement: If the values of a function belonging to the span of the reference shape functions are known in the evaluation nodes, then this function is uniquely determined. This entails that the number of evaluation nodes must be at least as big as the number of reference shape functions.
For triangular Lagrangian finite elements of degree p the evaluation nodes, usually called "interpolation nodes" in this context, can be chosen as \(\left(\frac{j}{p},\frac{k}{p}\right),\; 0\leq j,k \leq p, j+k\leq p\).
For some finite element spaces the interpolation functional may be defined based on integrals over edges. In this case the evaluation nodes will be quadrature nodes for the approximate evaluation of these integrals.
The quadrature rule must be exact for the polynomials contained in the local finite element spaces.
Implemented in lf::fe::FePoint< SCALAR >, lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Quad< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, lf::uscalfe::FeLagrangeO3Quad< SCALAR >, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >.
Referenced by lf::fe::ScalarReferenceFiniteElement< SCALAR >::PrintInfo().
|
pure virtual |
Computation of the gradients of all reference shape functions in a number of points.
refcoords | coordinates of N points in the reference cell passed as columns of a matrix of size dim x N. |
NumRefShapeFunctions() x (dim * refcoords.cols())
where dim
is the dimension of the reference finite element. The gradients are returned in chunks of rows of this matrix.Concerning the numbering of local shape functions, please consult the documentation of lf::assemble::DofHandler.
There are three reference shape functions \(\hat{b}^1,\hat{b}^2,\hat{b}^3\) associated with the vertices of the reference triangle. Let us assume that the refcoords
argument is a 2x2 matrix \([\mathbf{x}_1\;\mathbf{x}_2]\), which corresponds to passing the coordinates of two points in the reference triangle. Then this method will return a 3x4
matrix:
\[ \begin{pmatrix} \grad^T\hat{b}^1(\mathbf{x}_1) & \grad^T\hat{b}^1(\mathbf{x}_2) \\ \grad^T\hat{b}^2(\mathbf{x}_1) & \grad^T\hat{b}^2(\mathbf{x}_2) \\ \grad^T\hat{b}^3(\mathbf{x}_1) & \grad^T\hat{b}^3(\mathbf{x}_2) \end{pmatrix} \]
Implemented in lf::fe::FePoint< SCALAR >, lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >, lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Quad< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, and lf::uscalfe::FeLagrangeO3Quad< SCALAR >.
|
inlinevirtual |
Computes the linear combination of reference shape functions matching function values at evaluation nodes.
nodvals | row vector of function values at evaluation nodes The length of this vector must agree with NumEvaluationNodes(). |
If the evaluation nodes are interpolation nodes, that is, if the set of reference shape functions forms a cardinal basis with respect to these nodes, then we have NumEvaluationNodes() == NumRefShapeFunctions() and the linear mapping realized by NodalValuesToDofs() is the identity mapping.
If the vector of values at the evaluation nodes agrees with a vector of function values of a linear combination of reference shape functions at the evaluation nodes, then this method must return the very coefficients of the linear combination.
Reimplemented in lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >.
Definition at line 311 of file scalar_reference_finite_element.h.
References lf::fe::ScalarReferenceFiniteElement< SCALAR >::NumEvaluationNodes().
|
pure virtual |
Tell the number of evaluation (interpolation) nodes.
Implemented in lf::fe::FePoint< SCALAR >, lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Quad< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, lf::uscalfe::FeLagrangeO3Quad< SCALAR >, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >.
Referenced by lf::fe::ScalarReferenceFiniteElement< SCALAR >::NodalValuesToDofs(), and lf::fe::ScalarReferenceFiniteElement< SCALAR >::PrintInfo().
|
inlinevirtual |
Total number of reference shape functions associated with the reference cell.
Reimplemented in lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Quad< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, lf::uscalfe::FeLagrangeO3Quad< SCALAR >, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >.
Definition at line 128 of file scalar_reference_finite_element.h.
References lf::fe::ScalarReferenceFiniteElement< SCALAR >::Dimension(), lf::fe::ScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), lf::base::RefEl::NumSubEntities(), and lf::fe::ScalarReferenceFiniteElement< SCALAR >::RefEl().
Referenced by lf::fe::ScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), and lf::fe::ScalarReferenceFiniteElement< SCALAR >::PrintInfo().
|
inlinevirtual |
The number of interior reference shape functions for sub-entities of a particular co-dimension.
codim | co-dimension of the subentity |
Reimplemented in lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Quad< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, lf::uscalfe::FeLagrangeO3Quad< SCALAR >, lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >, and lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >.
Definition at line 152 of file scalar_reference_finite_element.h.
|
pure virtual |
The number of interior reference shape functions for every sub-entity.
codim | do-dimension of the subentity |
subidx | local index of the sub-entity |
Implemented in lf::fe::FePoint< SCALAR >, lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, lf::uscalfe::FeLagrangeO3Quad< SCALAR >, lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, and lf::uscalfe::FeLagrangeO2Quad< SCALAR >.
|
protecteddefault |
|
protecteddefaultnoexcept |
|
pure virtual |
Tells the type of reference cell underlying the parametric finite element.
Implemented in lf::fe::FePoint< SCALAR >, lf::fe::FeHierarchicSegment< SCALAR >, lf::fe::FeHierarchicTria< SCALAR >, lf::fe::FeHierarchicQuad< SCALAR >, lf::fe::test_utils::ComplexScalarReferenceFiniteElement< SCALAR >, lf::uscalfe::FeLagrangeO1Tria< SCALAR >, lf::uscalfe::FeLagrangeO1Quad< SCALAR >, lf::uscalfe::FeLagrangeO1Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Tria< SCALAR >, lf::uscalfe::FeLagrangeO2Segment< SCALAR >, lf::uscalfe::FeLagrangeO2Quad< SCALAR >, lf::uscalfe::FeLagrangeO3Tria< SCALAR >, lf::uscalfe::FeLagrangeO3Segment< SCALAR >, lf::uscalfe::FeLagrangeO3Quad< SCALAR >, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >.
Referenced by lf::fe::ScalarReferenceFiniteElement< SCALAR >::Dimension(), and lf::fe::ScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions().
|
related |
Stream output operator: just calls the ScalarReferenceFiniteElement::print() method.
Definition at line 353 of file scalar_reference_finite_element.h.
|
related |
Print information about a ScalarReferenceFiniteElement to the given stream.
o | stream to which output is to be sent |
srfe | The ScalarReferenceFiniteElement that should be printed. |
ctrl | controls the level of detail that is printed (see below) |
Definition at line 337 of file scalar_reference_finite_element.h.
References lf::fe::ScalarReferenceFiniteElement< SCALAR >::Degree(), lf::fe::ScalarReferenceFiniteElement< SCALAR >::EvaluationNodes(), lf::fe::ScalarReferenceFiniteElement< SCALAR >::NumEvaluationNodes(), and lf::fe::ScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions().