LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Helper class which stores a ScalarReferenceFiniteElement with precomputed values at the nodes of a quadrature rule. More...
#include <lf/uscalfe/uscalfe.h>
Public Member Functions | |
PrecomputedScalarReferenceFiniteElement ()=default | |
Default constructor which does not initialize this class at all (invalid state). If any method is called upon it, an error is thrown. | |
PrecomputedScalarReferenceFiniteElement (const PrecomputedScalarReferenceFiniteElement &)=delete | |
PrecomputedScalarReferenceFiniteElement (PrecomputedScalarReferenceFiniteElement &&) noexcept=default | |
PrecomputedScalarReferenceFiniteElement & | operator= (const PrecomputedScalarReferenceFiniteElement &)=delete |
PrecomputedScalarReferenceFiniteElement & | operator= (PrecomputedScalarReferenceFiniteElement &&) noexcept=default |
PrecomputedScalarReferenceFiniteElement (lf::fe::ScalarReferenceFiniteElement< SCALAR > const *fe, quad::QuadRule qr) | |
Main constructor performing precomputations. | |
bool | isInitialized () const |
Tells initialization status of object. | |
base::RefEl | RefEl () const override |
Tells the type of reference cell underlying the parametric finite element. | |
unsigned | Degree () const override |
Request the maximal polynomial degree of the basis functions in this finite element. | |
size_type | NumRefShapeFunctions () const override |
Total number of reference shape functions associated with the reference cell. | |
size_type | NumRefShapeFunctions (dim_t codim) const override |
The number of interior reference shape functions for sub-entities of a particular co-dimension. | |
size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const override |
The number of interior reference shape functions for every sub-entity. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | EvalReferenceShapeFunctions (const Eigen::MatrixXd &local) const override |
Evaluation of all reference shape functions in a number of points. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | GradientsReferenceShapeFunctions (const Eigen::MatrixXd &local) const override |
Computation of the gradients of all reference shape functions in a number of points. | |
Eigen::MatrixXd | EvaluationNodes () const override |
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom, nodal interpolation in the simplest case. | |
size_type | NumEvaluationNodes () const override |
Tell the number of evaluation (interpolation) nodes. | |
Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > | NodalValuesToDofs (const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > &nodvals) const override |
Computes the linear combination of reference shape functions matching function values at evaluation nodes. | |
const quad::QuadRule & | Qr () const |
Return the Quadrature rule at which the shape functions (and their gradients) have been precomputed. | |
const Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > & | PrecompReferenceShapeFunctions () const |
Value of EvalReferenceShapeFunctions(Qr().Points()) | |
const Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > & | PrecompGradientsReferenceShapeFunctions () const |
Value of EvalGradientsReferenceShapeFunctions(Qr().Points()) | |
~PrecomputedScalarReferenceFiniteElement () override=default | |
![]() | |
virtual | ~ScalarReferenceFiniteElement ()=default |
dim_t | Dimension () const |
Returns the spatial dimension of the reference cell. | |
Private Attributes | |
fe::ScalarReferenceFiniteElement< SCALAR > const * | fe_ = nullptr |
quad::QuadRule | qr_ |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | shap_fun_ |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | grad_shape_fun_ |
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. | |
Helper class which stores a ScalarReferenceFiniteElement with precomputed values at the nodes of a quadrature rule.
SCALAR | The scalar type of the shape functions, e.g. double |
This class does essentially three things:
Precomputing entity-independent quantities boost efficiency in the context of parametric finite element methods provided that a uniform quadrature rule is used for local computations on all mesh entities of the same topological type.
Detailed explanations can be found in Paragraph 2.8.3.6, Paragraph 2.7.5.8.
Definition at line 44 of file precomputed_scalar_reference_finite_element.h.
|
default |
Default constructor which does not initialize this class at all (invalid state). If any method is called upon it, an error is thrown.
|
delete |
|
defaultnoexcept |
|
inline |
Main constructor performing precomputations.
fe | definition of reference finite element |
qr | quadrature rule (nodes and weights on reference element) |
Initialization of local data members.
Definition at line 72 of file precomputed_scalar_reference_finite_element.h.
|
overridedefault |
virtual destructor
|
inlineoverridevirtual |
Request the maximal polynomial degree of the basis functions in this finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 92 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
inlineoverridevirtual |
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} \]
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 114 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
inlineoverridevirtual |
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.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 124 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
inlineoverridevirtual |
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} \]
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 119 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
inline |
Tells initialization status of object.
An object is in an undefined state when built by the default constructor
Definition at line 85 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval().
|
inlineoverridevirtual |
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 from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 133 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
inlineoverridevirtual |
Tell the number of evaluation (interpolation) nodes.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 128 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
inlineoverridevirtual |
Total number of reference shape functions associated with the reference cell.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 97 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval().
|
inlineoverridevirtual |
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 >.
Definition at line 102 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
inlineoverridevirtual |
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 >.
Definition at line 107 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
delete |
|
defaultnoexcept |
|
inline |
Value of EvalGradientsReferenceShapeFunctions(Qr().Points())
See fe::ScalarReferenceFiniteElement::GradientsReferenceShapeFunctions() for the packed format in which the gradients are returned.
Definition at line 166 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::grad_shape_fun_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval().
|
inline |
Value of EvalReferenceShapeFunctions(Qr().Points())
Definition at line 152 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::shap_fun_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval().
|
inline |
Return the Quadrature rule at which the shape functions (and their gradients) have been precomputed.
Definition at line 143 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_, and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::qr_.
Referenced by lf::uscalfe::ReactionDiffusionElementMatrixProvider< SCALAR, DIFF_COEFF, REACTION_COEFF >::Eval().
|
inlineoverridevirtual |
Tells the type of reference cell underlying the parametric finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 87 of file precomputed_scalar_reference_finite_element.h.
References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.
|
private |
The underlying scalar-valued parametric finite element
Definition at line 176 of file precomputed_scalar_reference_finite_element.h.
Referenced by lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::Degree(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::EvalReferenceShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::EvaluationNodes(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::GradientsReferenceShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::isInitialized(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NodalValuesToDofs(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumEvaluationNodes(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompGradientsReferenceShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompReferenceShapeFunctions(), lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::Qr(), and lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::RefEl().
|
private |
Holds gradients of reference shape functions at quadrature nodes
Definition at line 183 of file precomputed_scalar_reference_finite_element.h.
Referenced by lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompGradientsReferenceShapeFunctions().
|
private |
Uniform parametric quadrature rule for the associated type of reference element
Definition at line 179 of file precomputed_scalar_reference_finite_element.h.
Referenced by lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::Qr().
|
private |
Holds values of reference shape functions at reference quadrature nodes
Definition at line 181 of file precomputed_scalar_reference_finite_element.h.
Referenced by lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompReferenceShapeFunctions().