LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
lf::fe::HierarchicScalarFESpace< SCALAR > Class Template Reference

Finite Element Space that supports arbitrary, local degrees. More...

#include <lf/fe/fe.h>

Inheritance diagram for lf::fe::HierarchicScalarFESpace< SCALAR >:
lf::fe::ScalarFESpace< SCALAR >

Public Types

using Scalar = SCALAR
 
- Public Types inherited from lf::fe::ScalarFESpace< SCALAR >
using Scalar = SCALAR
 

Public Member Functions

 HierarchicScalarFESpace ()=delete
 
 HierarchicScalarFESpace (const HierarchicScalarFESpace &)=delete
 
 HierarchicScalarFESpace (HierarchicScalarFESpace &&) noexcept=default
 
HierarchicScalarFESpaceoperator= (const HierarchicScalarFESpace &)=delete
 
HierarchicScalarFESpaceoperator= (HierarchicScalarFESpace &&) noexcept=default
 
 HierarchicScalarFESpace (const std::shared_ptr< const lf::mesh::Mesh > &mesh_p, unsigned degree)
 Construct a new Hierarchic FESpace with uniform polynomial degree. More...
 
template<class F , class = std::enable_if_t< std::is_invocable_v<F, const mesh::Entity &>>>
 HierarchicScalarFESpace (const std::shared_ptr< const lf::mesh::Mesh > &mesh_p, F &&degree_functor)
 Construct a new Hierarchic FESpace with (possibly) varying polynomial degrees. More...
 
std::shared_ptr< const lf::mesh::MeshMesh () const override
 access to underlying mesh More...
 
const lf::assemble::DofHandlerLocGlobMap () const override
 access to associated local-to-global map More...
 
ScalarReferenceFiniteElement< SCALAR > const * ShapeFunctionLayout (const lf::mesh::Entity &entity) const override
 access to shape function layout for cells More...
 
size_type NumRefShapeFunctions (const lf::mesh::Entity &entity) const override
 number of interior shape functions associated to entities of various types More...
 
 ~HierarchicScalarFESpace () override=default
 
- Public Member Functions inherited from lf::fe::ScalarFESpace< SCALAR >
virtual std::shared_ptr< const lf::mesh::MeshMesh () const =0
 acess to underlying mesh More...
 
virtual const lf::assemble::DofHandlerLocGlobMap () const =0
 access to associated local-to-global map More...
 
virtual ScalarReferenceFiniteElement< SCALAR > const * ShapeFunctionLayout (const lf::mesh::Entity &entity) const =0
 access to shape function layout for mesh entities More...
 
virtual size_type NumRefShapeFunctions (const lf::mesh::Entity &entity) const =0
 number of interior shape functions associated to a particular mesh entity. More...
 
virtual ~ScalarFESpace ()=default
 No special destructor. More...
 

Private Member Functions

template<class F >
void Init (F &&degree_functor)
 

Private Attributes

std::shared_ptr< const lf::mesh::Meshmesh_p_
 
quad::QuadRuleCache qr_cache_
 
lf::mesh::utils::AllCodimMeshDataSet< std::variant< std::monostate, FePoint< SCALAR >, FeHierarchicSegment< SCALAR >, FeHierarchicTria< SCALAR >, FeHierarchicQuad< SCALAR > > > ref_el_
 
lf::assemble::DynamicFEDofHandler dofh_
 

Additional Inherited Members

- Protected Member Functions inherited from lf::fe::ScalarFESpace< SCALAR >
 ScalarFESpace ()=default
 default constructor, needed by std::vector More...
 
 ScalarFESpace (const ScalarFESpace &)=default
 
 ScalarFESpace (ScalarFESpace &&) noexcept=default
 
ScalarFESpaceoperator= (const ScalarFESpace &)=default
 
ScalarFESpaceoperator= (ScalarFESpace &&) noexcept=default
 

Detailed Description

template<typename SCALAR>
class lf::fe::HierarchicScalarFESpace< SCALAR >

Finite Element Space that supports arbitrary, local degrees.

Template Parameters
SCALARunderlying scalar type, usually either double or complex<double>

This FE Space contains hierarchic Finite Elements, meaning that the function spaces of lower polynomial degree are contained in the higher order function spaces.

The polynomial degree can vary from entity to entity, i.e. local \(p\)-refinement is supported.

Example usage

The following code snippet computes the solution of the BVP

\begin{align} - \Delta u &= 1 && \text{on }\Omega := [0,1]^2 \\ u &= 0 && \text{on }\partial \Omega \end{align}

// Get a mesh
// Create HierarchicalFESpace with uniform degree 5.
const unsigned degree = 5;
const auto fe_space =
std::make_shared<lf::fe::HierarchicScalarFESpace<double>>(mesh_p, degree);
// define diffusion coefficient
// define rhs load
// Assemble the system matrix and right hand side
Eigen::VectorXd rhs = Eigen::VectorXd::Zero(fe_space->LocGlobMap().NumDofs());
const assemble::DofHandler &dofh = fe_space->LocGlobMap();
assemble::COOMatrix<double> A_COO(dofh.NumDofs(), dofh.NumDofs());
DiffusionElementMatrixProvider element_matrix_provider(fe_space, mf_alpha);
AssembleMatrixLocally(0, dofh, dofh, element_matrix_provider, A_COO);
ScalarLoadElementVectorProvider element_vector_provider(fe_space, mf_load);
AssembleVectorLocally(0, dofh, element_vector_provider, rhs);
// Enforce zero dirichlet boundary conditions
const auto boundary = lf::mesh::utils::flagEntitiesOnBoundary(mesh_p);
const auto selector = [&](unsigned int idx) -> std::pair<bool, double> {
const auto &entity = dofh.Entity(idx);
return {entity.Codim() > 0 && boundary(entity), 0};
};
FixFlaggedSolutionComponents(selector, A_COO, rhs);
// Solve the LSE using the cholesky decomposition
Eigen::SparseMatrix<double> A = A_COO.makeSparse();
Eigen::SimplicialLDLT<Eigen::SparseMatrix<double>> solver(A);
const Eigen::VectorXd solution = solver.solve(rhs);
// visualize the solution:
io::VtkWriter vtk(mesh_p, "solution.vtk", 0, 5);
vtk.WritePointData("solution", MeshFunctionFE(fe_space, solution));
A MeshFunction which takes the same constant value on the whole mesh.
void AssembleMatrixLocally(dim_t codim, const DofHandler &dof_handler_trial, const DofHandler &dof_handler_test, ENTITY_MATRIX_PROVIDER &entity_matrix_provider, TMPMATRIX &matrix)
Assembly function for standard assembly of finite element matrices.
Definition: assembler.h:113
void AssembleVectorLocally(dim_t codim, const DofHandler &dof_handler, ENTITY_VECTOR_PROVIDER &entity_vector_provider, VECTOR &resultvector)
entity-local assembly of (right-hand-side) vectors from element vectors
Definition: assembler.h:297
void FixFlaggedSolutionComponents(SELECTOR &&selectvals, COOMatrix< SCALAR > &A, RHSVECTOR &b)
enforce prescribed solution components
Definition: fix_dof.h:87
ScalarLoadElementVectorProvider(PTR fe_space, MESH_FUNCTION mf) -> ScalarLoadElementVectorProvider< typename PTR::element_type::Scalar, MESH_FUNCTION >
MeshFunctionFE(std::shared_ptr< T >, const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &) -> MeshFunctionFE< typename T::Scalar, SCALAR_COEFF >
DiffusionElementMatrixProvider(PTR fe_space, DIFF_COEFF alpha) -> DiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF >
std::shared_ptr< lf::mesh::Mesh > GenerateHybrid2DTestMesh(int selector, double scale)
Generates a simple 2D hybrid test mesh.
Definition: test_meshes.cc:14
CodimMeshDataSet< bool > flagEntitiesOnBoundary(const std::shared_ptr< const Mesh > &mesh_p, lf::base::dim_t codim)
flag entities of a specific co-dimension located on the boundary

Definition at line 43 of file hierarchic_scalar_fe_space.h.

Member Typedef Documentation

◆ Scalar

template<typename SCALAR >
using lf::fe::HierarchicScalarFESpace< SCALAR >::Scalar = SCALAR

Definition at line 45 of file hierarchic_scalar_fe_space.h.

Constructor & Destructor Documentation

◆ HierarchicScalarFESpace() [1/5]

template<typename SCALAR >
lf::fe::HierarchicScalarFESpace< SCALAR >::HierarchicScalarFESpace ( )
delete

◆ HierarchicScalarFESpace() [2/5]

template<typename SCALAR >
lf::fe::HierarchicScalarFESpace< SCALAR >::HierarchicScalarFESpace ( const HierarchicScalarFESpace< SCALAR > &  )
delete

◆ HierarchicScalarFESpace() [3/5]

template<typename SCALAR >
lf::fe::HierarchicScalarFESpace< SCALAR >::HierarchicScalarFESpace ( HierarchicScalarFESpace< SCALAR > &&  )
defaultnoexcept

◆ HierarchicScalarFESpace() [4/5]

template<typename SCALAR >
lf::fe::HierarchicScalarFESpace< SCALAR >::HierarchicScalarFESpace ( const std::shared_ptr< const lf::mesh::Mesh > &  mesh_p,
unsigned  degree 
)
inlineexplicit

Construct a new Hierarchic FESpace with uniform polynomial degree.

Parameters
mesh_pA shared pointer to the underlying mesh (immutable)
degreeThe uniform polynomial degree.

Definition at line 59 of file hierarchic_scalar_fe_space.h.

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

◆ HierarchicScalarFESpace() [5/5]

template<typename SCALAR >
template<class F , class = std::enable_if_t< std::is_invocable_v<F, const mesh::Entity &>>>
lf::fe::HierarchicScalarFESpace< SCALAR >::HierarchicScalarFESpace ( const std::shared_ptr< const lf::mesh::Mesh > &  mesh_p,
F &&  degree_functor 
)
inlineexplicit

Construct a new Hierarchic FESpace with (possibly) varying polynomial degrees.

Template Parameters
Ftype of the degree_functor
Parameters
mesh_pA shared pointer to the underlying mesh (immutable)
degree_functorA function object that assigns a polynomial degree to every entity in the mesh. See below for more info.

Degree Functor

The degree_functor must overload the call operator as follows:

unsigned degree_functor(const lf::mesh::Entity& e)
Interface class representing a topological entity in a cellular complex
Definition: entity.h:39

and should return the polynomial degree for the respective entity. The degree functor will be called for all entities (i.e. edges, triangles, quadrilaterals) in the mesh except for points.

Definition at line 88 of file hierarchic_scalar_fe_space.h.

◆ ~HierarchicScalarFESpace()

template<typename SCALAR >
lf::fe::HierarchicScalarFESpace< SCALAR >::~HierarchicScalarFESpace ( )
overridedefault

Member Function Documentation

◆ Init()

template<typename SCALAR >
template<class F >
void lf::fe::HierarchicScalarFESpace< SCALAR >::Init ( F &&  degree_functor)
inlineprivate

◆ LocGlobMap()

template<typename SCALAR >
const lf::assemble::DofHandler & lf::fe::HierarchicScalarFESpace< SCALAR >::LocGlobMap ( ) const
inlineoverridevirtual

access to associated local-to-global map

Returns
a reference to the lf::assemble::DofHandler object (immutable)

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

Definition at line 124 of file hierarchic_scalar_fe_space.h.

◆ Mesh()

template<typename SCALAR >
std::shared_ptr< const lf::mesh::Mesh > lf::fe::HierarchicScalarFESpace< SCALAR >::Mesh ( ) const
inlineoverridevirtual

access to underlying mesh

Returns
a shared pointer to the mesh

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

Definition at line 116 of file hierarchic_scalar_fe_space.h.

◆ NumRefShapeFunctions()

template<typename SCALAR >
size_type lf::fe::HierarchicScalarFESpace< SCALAR >::NumRefShapeFunctions ( const lf::mesh::Entity entity) const
inlineoverridevirtual

number of interior shape functions associated to entities of various types

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

Definition at line 153 of file hierarchic_scalar_fe_space.h.

◆ operator=() [1/2]

template<typename SCALAR >
HierarchicScalarFESpace & lf::fe::HierarchicScalarFESpace< SCALAR >::operator= ( const HierarchicScalarFESpace< SCALAR > &  )
delete

◆ operator=() [2/2]

template<typename SCALAR >
HierarchicScalarFESpace & lf::fe::HierarchicScalarFESpace< SCALAR >::operator= ( HierarchicScalarFESpace< SCALAR > &&  )
defaultnoexcept

◆ ShapeFunctionLayout()

template<typename SCALAR >
ScalarReferenceFiniteElement< SCALAR > const * lf::fe::HierarchicScalarFESpace< SCALAR >::ShapeFunctionLayout ( const lf::mesh::Entity entity) const
inlineoverridevirtual

access to shape function layout for cells

access to shape function layout for mesh entities

Parameters
entityThe entity to get the reference element for
Warning
NULL pointers may be returned by this method in case a finite element specification was not given for a particular topological type of entity.
Note
The returned ShapeFunctionLayout pointer will remain for the entire lifetime of the owning ScalarFESpace.
See also
ScalarReferenceFiniteElement

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

Definition at line 132 of file hierarchic_scalar_fe_space.h.

Member Data Documentation

◆ dofh_

template<typename SCALAR >
lf::assemble::DynamicFEDofHandler lf::fe::HierarchicScalarFESpace< SCALAR >::dofh_
private

Definition at line 169 of file hierarchic_scalar_fe_space.h.

◆ mesh_p_

template<typename SCALAR >
std::shared_ptr<const lf::mesh::Mesh> lf::fe::HierarchicScalarFESpace< SCALAR >::mesh_p_
private

Definition at line 162 of file hierarchic_scalar_fe_space.h.

◆ qr_cache_

template<typename SCALAR >
quad::QuadRuleCache lf::fe::HierarchicScalarFESpace< SCALAR >::qr_cache_
private

Definition at line 163 of file hierarchic_scalar_fe_space.h.

◆ ref_el_

template<typename SCALAR >
lf::mesh::utils::AllCodimMeshDataSet< std::variant<std::monostate, FePoint<SCALAR>, FeHierarchicSegment<SCALAR>, FeHierarchicTria<SCALAR>, FeHierarchicQuad<SCALAR> > > lf::fe::HierarchicScalarFESpace< SCALAR >::ref_el_
private

Definition at line 168 of file hierarchic_scalar_fe_space.h.


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