LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
Quadratic Lagrangian finite element on a line segment. More...
#include <lf/uscalfe/uscalfe.h>
Public Member Functions | |
FeLagrangeO2Segment (const FeLagrangeO2Segment &)=default | |
FeLagrangeO2Segment (FeLagrangeO2Segment &&) noexcept=default | |
FeLagrangeO2Segment & | operator= (const FeLagrangeO2Segment &)=default |
FeLagrangeO2Segment & | operator= (FeLagrangeO2Segment &&) noexcept=default |
FeLagrangeO2Segment ()=default | |
~FeLagrangeO2Segment () 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 element rely on polynomials of degree 2. | |
size_type | NumRefShapeFunctions () const override |
Three local shape functions are associated with an edge. | |
size_type | NumRefShapeFunctions (dim_t codim) const override |
One shape function attached to each node and one interior shape function. | |
size_type | NumRefShapeFunctions (dim_t codim, sub_idx_t) const override |
One shape function attached to each node and one interior shape function. | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | EvalReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override |
Evaluation of shape functions on reference segment (unit interval) | |
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > | GradientsReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override |
Evaluation of derivatives of local shape function on reference interval. | |
Eigen::MatrixXd | EvaluationNodes () const override |
Evaluation nodes are the endpoints of the segment and its midpoint. | |
size_type | NumEvaluationNodes () const override |
Three evaluation nodes, same number as local shape functions. | |
![]() | |
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. | |
Quadratic Lagrangian finite element on a line segment.
This is a specialization of ScalarReferenceFiniteElement for quadratic Lagrangian finite elements on a segment.
The first two shape functions are associated with the vertices of the segment. The last one is an interior shape function. This complies with the ordering convention of LehreFEM++ for local shape functions, see Paragraph 2.7.4.11.
|
default |
|
defaultnoexcept |
|
default |
|
overridedefault |
|
inlineoverridevirtual |
Quadratic Lagrangian finite element rely on polynomials of degree 2.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Evaluation of shape functions on reference segment (unit interval)
The local shape functions are associated with a segment in the case of quadratic Lagrangian finite elements. In reference coordinates those are:
\[ \hat{b}^1(x) = 2(1-x)(\frac12-x)\;,\quad \hat{b}^2(x) = 2x(x-\frac12)\;,\quad \hat{b}^3(x) = 4(1-x)x\;. \]
This is a cardinal basis of the space of quadratic univariate polynomials with respect to the points \(x=0,\frac12,1\).
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Evaluation nodes are the endpoints of the segment and its midpoint.
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 >.
|
inlineoverridevirtual |
Evaluation of derivatives of local shape function on reference interval.
For the passed points return the derivative of the three local shape functions on the reference interval in the rows of a matrix.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
Three evaluation nodes, same number as local shape functions.
Tell the number of evaluation (interpolation) nodes.
Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 683 of file lagr_fe.h.
References lf::uscalfe::FeLagrangeO2Segment< SCALAR >::NumRefShapeFunctions().
|
inlineoverridevirtual |
Three local shape functions are associated with an edge.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
Definition at line 586 of file lagr_fe.h.
Referenced by lf::uscalfe::FeLagrangeO2Segment< SCALAR >::NumEvaluationNodes().
|
inlineoverridevirtual |
One shape function attached to each node and one interior shape function.
Hence (sub-)entities of any co-dimension possess exactly one local shape function.
Reimplemented from lf::fe::ScalarReferenceFiniteElement< SCALAR >.
|
inlineoverridevirtual |
One shape function attached to each node and one interior shape function.
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 573 of file lagr_fe.h.
References lf::base::RefEl::kSegment().