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

Quadratic Lagrangian finite element on a quadrilateral reference element. More...

#include <lf/uscalfe/uscalfe.h>

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

Public Member Functions

 FeLagrangeO2Quad (const FeLagrangeO2Quad &)=default
 
 FeLagrangeO2Quad (FeLagrangeO2Quad &&) noexcept=default
 
FeLagrangeO2Quadoperator= (const FeLagrangeO2Quad &)=default
 
FeLagrangeO2Quadoperator= (FeLagrangeO2Quad &&) noexcept=default
 
 FeLagrangeO2Quad ()=default
 
 ~FeLagrangeO2Quad () 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 elements sport polynomials of degree 2.
 
size_type NumRefShapeFunctions () const override
 Nine local shape functions are associated with a quadrilateral cell.
 
size_type NumRefShapeFunctions (dim_t codim) const override
 One shape function is attached to each node and edge of the quadrilateral. One interior shape function.
 
size_type NumRefShapeFunctions (dim_t codim, sub_idx_t) const override
 One shape function is attached to each node and each edge of the quadrilateral. One interior shape function.
 
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override
 Point evaluation of reference local shape functions in unit square.
 
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override
 Point evaluation of gradient of reference shape functions for quadratic Lagrangian finite elements on quadrilateral cells.
 
Eigen::MatrixXd EvaluationNodes () const override
 Evaluation nodes are the vertices, the midpoints of the edges and the center of the quadrilateral.
 
size_type NumEvaluationNodes () const override
 Nine evaluation nodes, same as the number of local shape functions.
 
- Public Member Functions inherited from lf::fe::ScalarReferenceFiniteElement< SCALAR >
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.
 

Static Private Attributes

static const FeLagrangeO2Segment< SCALAR > krsf_segment_
 
static const Eigen::Matrix< int, 9, 2 > ksegment_to_quad_mapping_
 

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<typename SCALAR>
class lf::uscalfe::FeLagrangeO2Quad< SCALAR >

Quadratic Lagrangian finite element on a quadrilateral reference element.

This is a specialization of ScalarReferenceFiniteElement. Refer to its documentation.

The shape functions are computed by a tensor product construction: The shape function associated with the interpolation node \( \mathbf{p} =(x_j,y_l)\) is computed as

\[ \hat{b}^{\mathbf{p}}(x_1,x_2)=\hat{b}^j(x_1)*\hat{b}^l(x_2) \]

Where \( \hat{b}^j\) and \( \hat{b}^l \) are the quadratic Lagrangian shape functions on the reference line segment associated to the interpolation nodes \( x_j \) and \( y_l \).

The first four shape functions are associated with the vertices, the next four with the edges of the quadrilateral. The last shape function is an interior shape function.

See also
ScalarReferenceFiniteElement

Definition at line 713 of file lagr_fe.h.

Constructor & Destructor Documentation

◆ FeLagrangeO2Quad() [1/3]

template<typename SCALAR >
lf::uscalfe::FeLagrangeO2Quad< SCALAR >::FeLagrangeO2Quad ( const FeLagrangeO2Quad< SCALAR > & )
default

◆ FeLagrangeO2Quad() [2/3]

template<typename SCALAR >
lf::uscalfe::FeLagrangeO2Quad< SCALAR >::FeLagrangeO2Quad ( FeLagrangeO2Quad< SCALAR > && )
defaultnoexcept

◆ FeLagrangeO2Quad() [3/3]

template<typename SCALAR >
lf::uscalfe::FeLagrangeO2Quad< SCALAR >::FeLagrangeO2Quad ( )
default

◆ ~FeLagrangeO2Quad()

template<typename SCALAR >
lf::uscalfe::FeLagrangeO2Quad< SCALAR >::~FeLagrangeO2Quad ( )
overridedefault

Member Function Documentation

◆ Degree()

template<typename SCALAR >
unsigned lf::uscalfe::FeLagrangeO2Quad< SCALAR >::Degree ( ) const
inlineoverridevirtual

Quadratic Lagrangian finite elements sport polynomials of degree 2.

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 730 of file lagr_fe.h.

◆ EvalReferenceShapeFunctions()

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

Point evaluation of reference local shape functions in unit square.

The nine local shape functions for quadratic Lagrangian finite elements on the unit square can be obtained by forming tensor product of the three local shape functions on the unit interval. This entails taking into account the numbering of the local shape functions: LSF 0-3 belong to the four corners, LSF 4-7 belong to the edges, LSF 8 belongs to the cell. This complies with the numbering convention for local shape functions as explained in Paragraph 2.7.4.11, Example 2.7.4.12.

Here the local shape functions are chosen as cardinal basis of the space of tensor product polynomials of degree 2 with respect to the following nine evaluation nodes: the four corners of the unit square, the four midpoints of its edges, and its center of gravity, see Figure 119.

See also
fe::ScalarReferenceFiniteElement::EvalReferenceShapeFunctions()

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 780 of file lagr_fe.h.

References lf::uscalfe::FeLagrangeO2Quad< SCALAR >::krsf_segment_, and lf::uscalfe::FeLagrangeO2Quad< SCALAR >::ksegment_to_quad_mapping_.

◆ EvaluationNodes()

template<typename SCALAR >
Eigen::MatrixXd lf::uscalfe::FeLagrangeO2Quad< SCALAR >::EvaluationNodes ( ) const
inlineoverridevirtual

Evaluation nodes are the vertices, the midpoints of the edges and the center of the quadrilateral.

Location and numbering of evaluation nodes:

Returns
a 2x9 matrix containing the refrence coordinates of the evaluation nodes in its columns

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 864 of file lagr_fe.h.

◆ GradientsReferenceShapeFunctions()

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

Point evaluation of gradient of reference shape functions for quadratic Lagrangian finite elements on quadrilateral cells.

Gradients are returned packed in the rows of a 9 x 2P matrix, where P is the number of evaluation points passed as a 2 x P matrix.

See also
FeLagrangeO2Quad::EvalReferenceShapeFunctions()
fe::ScalarReferenceFiniteElement::EvalReferenceShapeFunctions()

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 812 of file lagr_fe.h.

References lf::uscalfe::FeLagrangeO2Quad< SCALAR >::krsf_segment_, and lf::uscalfe::FeLagrangeO2Quad< SCALAR >::ksegment_to_quad_mapping_.

◆ NumEvaluationNodes()

template<typename SCALAR >
size_type lf::uscalfe::FeLagrangeO2Quad< SCALAR >::NumEvaluationNodes ( ) const
inlineoverridevirtual

Nine evaluation nodes, same as the number of local shape functions.

See also
FeLagrangeO2Quad::EvaluationNodes()

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 878 of file lagr_fe.h.

References lf::uscalfe::FeLagrangeO2Quad< SCALAR >::NumRefShapeFunctions().

◆ NumRefShapeFunctions() [1/3]

template<typename SCALAR >
size_type lf::uscalfe::FeLagrangeO2Quad< SCALAR >::NumRefShapeFunctions ( ) const
inlineoverridevirtual

Nine local shape functions are associated with a quadrilateral cell.

Refer to Example 2.6.2.7.

See also
fe::ScalarReferenceFiniteElement::NumRefShapeFunctions()

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

Definition at line 738 of file lagr_fe.h.

Referenced by lf::uscalfe::FeLagrangeO2Quad< SCALAR >::NumEvaluationNodes().

◆ NumRefShapeFunctions() [2/3]

template<typename SCALAR >
size_type lf::uscalfe::FeLagrangeO2Quad< SCALAR >::NumRefShapeFunctions ( dim_t codim) const
inlineoverridevirtual

One shape function is attached to each node and edge of the quadrilateral. One interior shape function.

See also
fe::ScalarReferenceFiniteElement::NumRefShapeFunctions(dim_t codim)

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

Definition at line 745 of file lagr_fe.h.

◆ NumRefShapeFunctions() [3/3]

template<typename SCALAR >
size_type lf::uscalfe::FeLagrangeO2Quad< SCALAR >::NumRefShapeFunctions ( dim_t codim,
sub_idx_t  ) const
inlineoverridevirtual

One shape function is attached to each node and each edge of the quadrilateral. One interior shape function.

See also
fe::ScalarReferenceFiniteElement::NumRefShapeFunctions(dim_t, sub_idx_t)

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 756 of file lagr_fe.h.

◆ operator=() [1/2]

template<typename SCALAR >
FeLagrangeO2Quad & lf::uscalfe::FeLagrangeO2Quad< SCALAR >::operator= ( const FeLagrangeO2Quad< SCALAR > & )
default

◆ operator=() [2/2]

template<typename SCALAR >
FeLagrangeO2Quad & lf::uscalfe::FeLagrangeO2Quad< SCALAR >::operator= ( FeLagrangeO2Quad< SCALAR > && )
defaultnoexcept

◆ RefEl()

template<typename SCALAR >
base::RefEl lf::uscalfe::FeLagrangeO2Quad< SCALAR >::RefEl ( ) const
inlineoverridevirtual

Tells the type of reference cell underlying the parametric finite element.

Implements lf::fe::ScalarReferenceFiniteElement< SCALAR >.

Definition at line 724 of file lagr_fe.h.

References lf::base::RefEl::kQuad().

Member Data Documentation

◆ krsf_segment_

template<typename SCALAR >
const FeLagrangeO2Segment< SCALAR > lf::uscalfe::FeLagrangeO2Quad< SCALAR >::krsf_segment_
staticprivate
Initial value:
=
FeLagrangeO2Segment<SCALAR>()

description of the shape functions on the reference line segment \( [0,1] \) that are used in the tensor product construction of the shape functions

Definition at line 886 of file lagr_fe.h.

Referenced by lf::uscalfe::FeLagrangeO2Quad< SCALAR >::EvalReferenceShapeFunctions(), and lf::uscalfe::FeLagrangeO2Quad< SCALAR >::GradientsReferenceShapeFunctions().

◆ ksegment_to_quad_mapping_

template<typename SCALAR >
const Eigen::Matrix< int, 9, 2 > lf::uscalfe::FeLagrangeO2Quad< SCALAR >::ksegment_to_quad_mapping_
staticprivate
Initial value:
=
(Eigen::Matrix<int,9,2>() << 0, 0,
1, 0,
1, 1,
0, 1,
2, 0,
1, 2,
2, 1,
0, 2,
2, 2).finished()

The shape functions are computed by a tensor product construction: The shape function associated with the interpolation node \( \mathbf{p} =(x_j,y_l)\) is computed as

\[ \hat{b}^{\mathbf{p}}(x_1,x_2)=\hat{b}^j(x_1)*\hat{b}^l(x_2) \]

Where \( \hat{b}^j\) and \( \hat{b}^l \) are the quadratic Lagrangian shape functions on the reference line segment associated to the interpolation nodes \( x_j \) and \( y_l \).

This matrix of size NumRefShapeFunctions x 2 contains in each row the indices of the "segment" shape functions, that define the reference shape function via the above formula. if \(p\) was the \(k\)-th interpolation node, the \(k\)-th row would contain the indices \(j\) and \(l\).

Definition at line 905 of file lagr_fe.h.

Referenced by lf::uscalfe::FeLagrangeO2Quad< SCALAR >::EvalReferenceShapeFunctions(), and lf::uscalfe::FeLagrangeO2Quad< SCALAR >::GradientsReferenceShapeFunctions().


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