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

Cubic Lagrangian finite elment on a triangular reference element. More...

#include <lf/uscalfe/uscalfe.h>

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

Public Member Functions

 FeLagrangeO3Tria (const FeLagrangeO3Tria &)=default
 
 FeLagrangeO3Tria (FeLagrangeO3Tria &&) noexcept=default
 
FeLagrangeO3Triaoperator= (const FeLagrangeO3Tria &)=default
 
FeLagrangeO3Triaoperator= (FeLagrangeO3Tria &&) noexcept=default
 
 FeLagrangeO3Tria ()=default
 
 ~FeLagrangeO3Tria () 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 rely on polynomials of degree 3.
 
size_type NumRefShapeFunctions () const override
 Ten local shape functions are associated with a triangular cell in the case of cubic Lagrangian finite elements.
 
size_type NumRefShapeFunctions (dim_t codim) const override
 One shape function attached to each node and two to each edge of the triangle. There is one interior shape function.
 
size_type NumRefShapeFunctions (dim_t codim, sub_idx_t subidx) const override
 One shape function attached to each node and two to each edge of the triangle. There is one interior shape function.
 
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override
 Point evaluation of reference shape functions.
 
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > GradientsReferenceShapeFunctions (const Eigen::MatrixXd &refcoords) const override
 Point evaluations of gradients of reference shape functions.
 
Eigen::MatrixXd EvaluationNodes () const override
 Evaluation nodes are the three vertices of the triangle, two equispaced interio nodes for each edge, a the barycentre.
 
size_type NumEvaluationNodes () const override
 Ten evaluation nodes.
 
- 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.
 

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::FeLagrangeO3Tria< SCALAR >

Cubic Lagrangian finite elment on a triangular reference element.

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

The 10 reference shape functions are combinations of the barycentric coordinate functions on the reference triangle.

The first three shape functions are associated with vertices, the next six with edges of the triangle. The last shape function is an interior shape function, see Example 2.6.1.7 for the location of interpolation nodes.

See also
ScalarReferenceFiniteElement

Definition at line 945 of file lagr_fe.h.

Constructor & Destructor Documentation

◆ FeLagrangeO3Tria() [1/3]

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

◆ FeLagrangeO3Tria() [2/3]

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

◆ FeLagrangeO3Tria() [3/3]

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

◆ ~FeLagrangeO3Tria()

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

Member Function Documentation

◆ Degree()

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

Quadratic Lagrangian finite elements rely on polynomials of degree 3.

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

Definition at line 962 of file lagr_fe.h.

◆ EvalReferenceShapeFunctions()

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

Point evaluation of reference shape functions.

The reference shape function are cardinal basis functions of the space of cubic polynomials with respect to the ten interpolation nodes supplied by the lattive drawn in Figure 111.

See also
fe::ScalarReferenceFiniteElement::EvalReferenceShapeFunctions()

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

Definition at line 1010 of file lagr_fe.h.

◆ EvaluationNodes()

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

Evaluation nodes are the three vertices of the triangle, two equispaced interio nodes for each edge, a the barycentre.

The reference shape functions as implemented by this class are a cardinal basis of the space of cubic two-variate polynomials with respect to these evaluation nodes.

The numbering of evluation nodes is fixed by the numbering of the local shape functions.

See also
fe::ScalarReferenceFiniteElement::EvaluationNodes()

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

Definition at line 1115 of file lagr_fe.h.

◆ GradientsReferenceShapeFunctions()

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

Point evaluations of gradients of reference shape functions.

Compute values of gradients of reference shape functions at specified nodes

See also
fe::ScalarReferenceFiniteElement::GradientsReferenceShapeFunctions()

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

Definition at line 1053 of file lagr_fe.h.

◆ NumEvaluationNodes()

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

Ten evaluation nodes.

Tell the number of evaluation (interpolation) nodes.

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

Definition at line 1137 of file lagr_fe.h.

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

◆ NumRefShapeFunctions() [1/3]

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

Ten local shape functions are associated with a triangular cell in the case of cubic Lagrangian finite elements.

See also
fe::ScalarReferenceFiniteElement::NumRefShapeFunctions()

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

Definition at line 969 of file lagr_fe.h.

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

◆ NumRefShapeFunctions() [2/3]

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

One shape function attached to each node and two to each edge of the triangle. There is one interior shape function.

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 979 of file lagr_fe.h.

◆ NumRefShapeFunctions() [3/3]

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

One shape function attached to each node and two to each edge of the triangle. There is one interior shape function.

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 992 of file lagr_fe.h.

◆ operator=() [1/2]

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

◆ operator=() [2/2]

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

◆ RefEl()

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

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

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

Definition at line 955 of file lagr_fe.h.

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


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