LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Member Functions | Private Attributes | List of all members
lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR > Class Template Reference

Helper class which stores a ScalarReferenceFiniteElement with precomputed values at the nodes of a quadrature rule. More...

#include <lf/uscalfe/uscalfe.h>

Inheritance diagram for lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >:
lf::fe::ScalarReferenceFiniteElement< SCALAR >

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
 
PrecomputedScalarReferenceFiniteElementoperator= (const PrecomputedScalarReferenceFiniteElement &)=delete
 
PrecomputedScalarReferenceFiniteElementoperator= (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::QuadRuleQr () 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
 
- Public Member Functions inherited from lf::fe::ScalarReferenceFiniteElement< SCALAR >
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

- Public Types inherited from lf::fe::ScalarReferenceFiniteElement< SCALAR >
using Scalar = SCALAR
 The underlying scalar type.
 
- Protected Member Functions inherited from lf::fe::ScalarReferenceFiniteElement< SCALAR >
 ScalarReferenceFiniteElement ()=default
 
 ScalarReferenceFiniteElement (const ScalarReferenceFiniteElement &)=default
 
 ScalarReferenceFiniteElement (ScalarReferenceFiniteElement &&) noexcept=default
 
ScalarReferenceFiniteElementoperator= (const ScalarReferenceFiniteElement &)=default
 
ScalarReferenceFiniteElementoperator= (ScalarReferenceFiniteElement &&) noexcept=default
 

Detailed Description

template<class SCALAR>
class lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >

Helper class which stores a ScalarReferenceFiniteElement with precomputed values at the nodes of a quadrature rule.

Template Parameters
SCALARThe scalar type of the shape functions, e.g. double

This class does essentially three things:

  1. It wraps any ScalarReferenceFiniteElement and forwards all calls to the wrapped instance.
  2. It provides access to a quad::QuadratureRule (passed in constructor)
  3. It provides additional member functions to access the precomputed values
  4. the shape functions/gradients at the nodes of the quadrature rule.

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.

Constructor & Destructor Documentation

◆ PrecomputedScalarReferenceFiniteElement() [1/4]

template<class SCALAR >
lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::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() [2/4]

template<class SCALAR >
lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecomputedScalarReferenceFiniteElement ( const PrecomputedScalarReferenceFiniteElement< SCALAR > & )
delete

◆ PrecomputedScalarReferenceFiniteElement() [3/4]

template<class SCALAR >
lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecomputedScalarReferenceFiniteElement ( PrecomputedScalarReferenceFiniteElement< SCALAR > && )
defaultnoexcept

◆ PrecomputedScalarReferenceFiniteElement() [4/4]

template<class SCALAR >
lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecomputedScalarReferenceFiniteElement ( lf::fe::ScalarReferenceFiniteElement< SCALAR > const * fe,
quad::QuadRule qr )
inline

Main constructor performing precomputations.

Parameters
fedefinition of reference finite element
qrquadrature rule (nodes and weights on reference element)

Initialization of local data members.

Definition at line 72 of file precomputed_scalar_reference_finite_element.h.

◆ ~PrecomputedScalarReferenceFiniteElement()

virtual destructor

Member Function Documentation

◆ Degree()

template<class SCALAR >
unsigned lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::Degree ( ) const
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_.

◆ EvalReferenceShapeFunctions()

template<class SCALAR >
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::EvalReferenceShapeFunctions ( const Eigen::MatrixXd & refcoords) const
inlineoverridevirtual

Evaluation of all reference shape functions in a number of points.

Parameters
refcoordscoordinates 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
Returns
An Eigen Matrix of size 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.

Note
shape functions are assumed to be real-valued.

Example: Triangular Linear Lagrangian finite elements

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_.

◆ EvaluationNodes()

template<class SCALAR >
Eigen::MatrixXd lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::EvaluationNodes ( ) const
inlineoverridevirtual

Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom, nodal interpolation in the simplest case.

Returns
A d x N matrix, where d is the dimension of the reference cell, and N the number of interpolation nodes. The columns of this matrix contain their reference coordinates.

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.

Unisolvence

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.

Note
It is not required that any vector a values at evaluation nodes can be produced by a suitable linear combination of reference shape functions. For instance, this will not be possible, if there are more evaluation points than reference shape functions. If both sets have the same size, however, the interpolation problem has a solution for any vector of values at the evluation points.

Example: Principal lattice

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\).

Moment-based degrees of freedom

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_.

◆ GradientsReferenceShapeFunctions()

template<class SCALAR >
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::GradientsReferenceShapeFunctions ( const Eigen::MatrixXd & refcoords) const
inlineoverridevirtual

Computation of the gradients of all reference shape functions in a number of points.

Parameters
refcoordscoordinates of N points in the reference cell passed as columns of a matrix of size dim x N.
Returns
An Eigen Matrix of size 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.

Example: Triangular Linear Lagrangian finite elements

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_.

◆ isInitialized()

template<class SCALAR >
bool lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::isInitialized ( ) const
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().

◆ NodalValuesToDofs()

template<class SCALAR >
Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NodalValuesToDofs ( const Eigen::Matrix< SCALAR, 1, Eigen::Dynamic > & nodvals) const
inlineoverridevirtual

Computes the linear combination of reference shape functions matching function values at evaluation nodes.

Parameters
nodvalsrow vector of function values at evaluation nodes The length of this vector must agree with NumEvaluationNodes().
Returns
The coefficients of the local "nodal interpolant" with respect to the reference shape functions. This is a row vector of length NumRefShapeFunctions().

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.

Note
default implementation is the identity mapping

Requirement: reproduction of local finite element functions

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_.

◆ NumEvaluationNodes()

template<class SCALAR >
size_type lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumEvaluationNodes ( ) const
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_.

◆ NumRefShapeFunctions() [1/3]

template<class SCALAR >
size_type lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions ( ) const
inlineoverridevirtual

Total number of reference shape functions associated with the reference cell.

Note
the total number of shape functions is the sum of the number of interior shape functions for all sub-entities and the entity itself.

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().

◆ NumRefShapeFunctions() [2/3]

template<class SCALAR >
size_type lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions ( dim_t codim) const
inlineoverridevirtual

The number of interior reference shape functions for sub-entities of a particular co-dimension.

Parameters
codimco-dimension of the subentity
Returns
number of interior reference shape function belonging to the specified sub-entity.
Note
this method will throw an exception when different numbers of reference shape functions are assigned to different sub-entities of the same co-dimension

Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 102 of file precomputed_scalar_reference_finite_element.h.

References lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_.

◆ NumRefShapeFunctions() [3/3]

template<class SCALAR >
size_type lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions ( dim_t codim,
sub_idx_t subidx ) const
inlineoverridevirtual

The number of interior reference shape functions for every sub-entity.

Parameters
codimdo-dimension of the subentity
subidxlocal index of the sub-entity
Returns
number of interior reference shape function belonging to the specified 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_.

◆ operator=() [1/2]

◆ operator=() [2/2]

◆ PrecompGradientsReferenceShapeFunctions()

template<class SCALAR >
const Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > & lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompGradientsReferenceShapeFunctions ( ) const
inline

◆ PrecompReferenceShapeFunctions()

template<class SCALAR >
const Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > & lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::PrecompReferenceShapeFunctions ( ) const
inline

◆ Qr()

template<class SCALAR >
const quad::QuadRule & lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::Qr ( ) const
inline

◆ RefEl()

template<class SCALAR >
base::RefEl lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::RefEl ( ) const
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_.

Member Data Documentation

◆ fe_

template<class SCALAR >
fe::ScalarReferenceFiniteElement<SCALAR> const* lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::fe_ = nullptr
private

◆ grad_shape_fun_

template<class SCALAR >
Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::grad_shape_fun_
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().

◆ qr_

template<class SCALAR >
quad::QuadRule lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::qr_
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().

◆ shap_fun_

template<class SCALAR >
Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> lf::uscalfe::PrecomputedScalarReferenceFiniteElement< SCALAR >::shap_fun_
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().


The documentation for this class was generated from the following file: