LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Decorator class around a ScalarReferenceFiniteElement to represent traces on the cell boundary. More...
#include <projects/dpg/trace_scalar_reference_finite_element.h>
Public Member Functions | |
TraceScalarReferenceFiniteElement ()=default | |
Default constructor, does not initialize this class (invalid state). More... | |
TraceScalarReferenceFiniteElement (const TraceScalarReferenceFiniteElement &)=delete | |
TraceScalarReferenceFiniteElement (TraceScalarReferenceFiniteElement &&) noexcept=default | |
TraceScalarReferenceFiniteElement & | operator= (const TraceScalarReferenceFiniteElement &)=delete |
TraceScalarReferenceFiniteElement & | operator= (TraceScalarReferenceFiniteElement &&) noexcept=default |
TraceScalarReferenceFiniteElement (std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > cfe) | |
bool | isInitialized () const |
Report the initialization status of the object. More... | |
lf::base::RefEl | RefEl () const override |
Tells the type of reference cell underlying the parametric finite element. More... | |
unsigned | Degree () const override |
Request the maximal polynomial degree of the basis functions in this finite element. More... | |
size_type | NumRefShapeFunctions () const override |
Total number of reference shape functions associated with the reference cell. More... | |
size_type | NumRefShapeFunctions (dim_t codim) const override |
The number of interior reference shape functions for sub-entities of a particular co-dimension. More... | |
size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const override |
The number of interior reference shape functions for every sub-entity. More... | |
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. More... | |
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. More... | |
Eigen::MatrixXd | EvaluationNodes () const override |
Returns reference coordinates of "evaluation nodes" for evaluation of parametric degrees of freedom, nodal interpolation in the simplest case. More... | |
lf::fe::size_type | NumEvaluationNodes () const override |
Tell the number of evaluation (interpolation) nodes. More... | |
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. More... | |
~TraceScalarReferenceFiniteElement () override=default | |
![]() | |
virtual | ~ScalarReferenceFiniteElement ()=default |
virtual base::RefEl | RefEl () const =0 |
Tells the type of reference cell underlying the parametric finite element. More... | |
virtual unsigned int | Degree () const =0 |
Request the maximal polynomial degree of the basis functions in this finite element. More... | |
dim_t | Dimension () const |
Returns the spatial dimension of the reference cell. More... | |
virtual size_type | NumRefShapeFunctions () const |
Total number of reference shape functions associated with the reference cell. More... | |
virtual size_type | NumRefShapeFunctions (dim_t codim) const |
The number of interior reference shape functions for sub-entities of a particular co-dimension. More... | |
virtual size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const =0 |
The number of interior reference shape functions for every sub-entity. More... | |
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. More... | |
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. More... | |
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. More... | |
virtual size_type | NumEvaluationNodes () const =0 |
Tell the number of evaluation (interpolation) nodes. More... | |
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. More... | |
Private Attributes | |
std::shared_ptr< const lf::fe::ScalarReferenceFiniteElement< SCALAR > > | cfe_ |
Additional Inherited Members | |
![]() | |
using | Scalar = SCALAR |
The underlying scalar type. More... | |
![]() | |
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. More... | |
template<typename SCALAR > | |
std::ostream & | operator<< (std::ostream &o, const ScalarReferenceFiniteElement< SCALAR > &fe_desc) |
Stream output operator: just calls the ScalarReferenceFiniteElement::print() method. More... | |
Decorator class around a ScalarReferenceFiniteElement to represent traces on the cell boundary.
SCALAR | The scalar type of the shape functions e.g.'double' |
The class decorates any lf::fe::ScalarReferenceFiniteElement and drops all its interior shape functions.
This class is used to represent finite elements discretizing functions in \( H^{1/2}(\mathcal S) \), the space of traces of functions in \( H^1(\Omega) \) on the mesh skeleton. Traces are evaluated/integrated only on/over the boundary of cells. If a cardinal basis is used in the underlying finite element, interior shape functions evalute to zero on any point of the cell boundary. To construct a basis corresponding to the traces they thus have to be dropped.
Definition at line 42 of file trace_scalar_reference_finite_element.h.
|
default |
Default constructor, does not initialize this class (invalid state).
If any method is called upon it, an error is thrown.
|
delete |
|
defaultnoexcept |
|
inlineexplicit |
Definition at line 61 of file trace_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 77 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized().
|
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 99 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized(), and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions().
|
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 113 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized(), and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions().
|
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 106 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized(), and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions().
|
inline |
Report the initialization status of the object.
Objects built by the default constructor are undefined
Definition at line 70 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_.
Referenced by projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::Degree(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::EvalReferenceShapeFunctions(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::EvaluationNodes(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::GradientsReferenceShapeFunctions(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NodalValuesToDofs(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumEvaluationNodes(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::RefEl().
|
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 123 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumEvaluationNodes(), and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions().
|
inlineoverridevirtual |
Tell the number of evaluation (interpolation) nodes.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 118 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized().
Referenced by projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NodalValuesToDofs().
|
inlineoverridevirtual |
Total number of reference shape functions associated with the reference cell.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 82 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized().
Referenced by projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::EvalReferenceShapeFunctions(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::EvaluationNodes(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::GradientsReferenceShapeFunctions(), and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NodalValuesToDofs().
|
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 87 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized().
|
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 92 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized().
|
delete |
|
defaultnoexcept |
|
inlineoverridevirtual |
Tells the type of reference cell underlying the parametric finite element.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 72 of file trace_scalar_reference_finite_element.h.
References projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::cfe_, and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized().
|
private |
The underlying scalar-valued paramteric finite element
Definition at line 147 of file trace_scalar_reference_finite_element.h.
Referenced by projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::Degree(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::EvalReferenceShapeFunctions(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::EvaluationNodes(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::GradientsReferenceShapeFunctions(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::isInitialized(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NodalValuesToDofs(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumEvaluationNodes(), projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::NumRefShapeFunctions(), and projects::dpg::TraceScalarReferenceFiniteElement< SCALAR >::RefEl().