LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Namespaces | Classes | Typedefs | Functions | Variables
lf::fe Namespace Reference

Collects data structures and algorithms designed for scalar finite element methods primarily meant for second-order elliptic boundary value problems. More...

Namespaces

namespace  test_utils
 Includes utilities to test classes from lf::fe.
 

Classes

class  DiffusionElementMatrixProvider
 Class for computing element matrices for general scalar-valued finite elements and homogeneous 2nd-order elliptic bilinear forms. More...
 
class  FeHierarchicQuad
 Hierarchic Finite Elements of arbitrary degree on quadrilaterals. More...
 
class  FeHierarchicSegment
 Hierarchic Finite Elements of arbitrary degree on segments. More...
 
class  FeHierarchicTria
 Hierarchic Finite Elements of arbitrary degree on triangles. More...
 
class  FePoint
 Linear finite element on a point. More...
 
class  HierarchicScalarFESpace
 Finite Element Space that supports arbitrary, local degrees. More...
 
class  MassEdgeMatrixProvider
 Quadrature-based computation of local mass matrix for an edge. More...
 
class  MassElementMatrixProvider
 Class for local quadrature based computation of element matrix for Lagrangian finite elements and a weighted \(L^2\) inner product. More...
 
class  MeshFunctionFE
 A MeshFunction representing an element from a ScalarUniformFESpace (e.g. solution of BVP) More...
 
class  MeshFunctionGradFE
 A MeshFunction representing the gradient of a function from a finite element space (e.g. gradient of a solution of BVP). More...
 
class  ScalarFESpace
 Space of scalar valued finite element functions on a hybrid 2D mesh More...
 
class  ScalarLoadEdgeVectorProvider
 Local edge contributions to element vector. More...
 
class  ScalarLoadElementVectorProvider
 Local computation of general element (load) vector for scalar finite elements; volume contributions only. More...
 
class  ScalarReferenceFiniteElement
 Interface class for parametric scalar valued finite elements. More...
 

Typedefs

using gdof_idx_t = lf::assemble::gdof_idx_t
 
using ldof_idx_t = lf::assemble::ldof_idx_t
 
using size_type = lf::assemble::size_type
 
using dim_t = lf::assemble::dim_t
 
using glb_idx_t = lf::assemble::glb_idx_t
 
using sub_idx_t = lf::base::sub_idx_t
 

Functions

template<class MF , class QR_SELECTOR , class ENTITY_PREDICATE = base::PredicateTrue, class = std::enable_if_t< std::is_invocable_v<QR_SELECTOR, const mesh::Entity &>>>
auto IntegrateMeshFunction (const lf::mesh::Mesh &mesh, const MF &mf, const QR_SELECTOR &qr_selector, const ENTITY_PREDICATE &ep=base::PredicateTrue{}, int codim=0) -> mesh::utils::MeshFunctionReturnType< MF >
 Integrate a MeshFunction over a mesh (with quadrature rules) More...
 
template<class MF , class ENTITY_PREDICATE = base::PredicateTrue>
auto IntegrateMeshFunction (const lf::mesh::Mesh &mesh, const MF &mf, int quad_degree, const ENTITY_PREDICATE &ep=base::PredicateTrue{}, int codim=0)
 Integrate a mesh function over a mesh using quadrature rules of uniform order. More...
 
template<typename SCALAR , typename MF , typename SELECTOR = base::PredicateTrue>
auto NodalProjection (const lf::fe::ScalarFESpace< SCALAR > &fe_space, MF &&u, SELECTOR &&pred=base::PredicateTrue{})
 Computes nodal projection of a mesh function and returns the finite element basis expansion coefficients of the result. More...
 
template<typename SCALAR , typename EDGESELECTOR , typename FUNCTION >
std::vector< std::pair< bool, SCALAR > > InitEssentialConditionFromFunction (const lf::fe::ScalarFESpace< SCALAR > &fes, EDGESELECTOR &&esscondflag, FUNCTION &&g)
 Initialization of flags/values for dofs of a Scalar finite element space whose values are imposed by a specified function. More...
 
double legendre (unsigned n, double x, double t)
 computes the n-th degree scaled Legendre Polynomial \( P_n(x;t) *\) More...
 
double ilegendre (unsigned n, double x, double t)
 computes the integral of the (n-1)-th degree scaled Legendre Polynomial More...
 
double ilegendre_dx (unsigned n, double x, double t)
 Computes \( \frac{\partial}{\partial x} L(x;t) \). More...
 
double ilegendre_dt (unsigned n, double x, double t)
 Computes \( \frac{\partial}{\partial t} L(x;t) \). More...
 
double legendre_dx (unsigned n, double x, double t)
 Computes the derivative of the n-th degree scaled Legendre polynomial. More...
 
double jacobi (unsigned n, double alpha, double beta, double x)
 Computes the n-th degree shifted Jacobi polynomial. More...
 
double jacobi (unsigned n, double alpha, double x)
 Computes the n-th degree shifted Jacobi polynomial for \( \beta = 0 \). More...
 
double ijacobi (unsigned n, double alpha, double x)
 Evaluate the integral of the (n-1)-th degree Jacobi Polynomial for \( *\beta = 0 \). More...
 
double ijacobi_dx (unsigned n, double alpha, double x)
 Computes the derivative of the n-th integrated scaled Jacobi polynomial. More...
 
double jacobi_dx (unsigned n, double alpha, double x)
 Computes the derivative of the n-th degree Jacobi Polynomial for \( *\beta = 0 \). More...
 
template<class PTR , class DIFF_COEFF >
 DiffusionElementMatrixProvider (PTR fe_space, DIFF_COEFF alpha) -> DiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF >
 
template<class PTR , class REACTION_COEFF >
 MassElementMatrixProvider (PTR fe_space, REACTION_COEFF gamma) -> MassElementMatrixProvider< typename PTR::element_type::Scalar, REACTION_COEFF >
 
template<class PTR , class COEFF , class EDGESELECTOR = base::PredicateTrue>
 MassEdgeMatrixProvider (PTR, COEFF coeff, EDGESELECTOR edge_predicate=base::PredicateTrue{}) -> MassEdgeMatrixProvider< typename PTR::element_type::Scalar, COEFF, EDGESELECTOR >
 
template<class PTR , class MESH_FUNCTION >
 ScalarLoadElementVectorProvider (PTR fe_space, MESH_FUNCTION mf) -> ScalarLoadElementVectorProvider< typename PTR::element_type::Scalar, MESH_FUNCTION >
 
template<class PTR , class FUNCTOR , class EDGESELECTOR = base::PredicateTrue>
 ScalarLoadEdgeVectorProvider (PTR, FUNCTOR, EDGESELECTOR=base::PredicateTrue{}) -> ScalarLoadEdgeVectorProvider< typename PTR::element_type::Scalar, FUNCTOR, EDGESELECTOR >
 
template<class T , class SCALAR_COEFF >
 MeshFunctionFE (std::shared_ptr< T >, const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &) -> MeshFunctionFE< typename T::Scalar, SCALAR_COEFF >
 
template<class T , class SCALAR_COEFF >
 MeshFunctionGradFE (std::shared_ptr< T >, const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &) -> MeshFunctionGradFE< typename T::Scalar, SCALAR_COEFF >
 
template<typename SCALAR_COEFF , typename FES_COARSE , typename FES_FINE >
Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > prolongate (const lf::refinement::MeshHierarchy &mh, std::shared_ptr< FES_COARSE > fespace_coarse, std::shared_ptr< FES_FINE > fespace_fine, const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &dofs_coarse, lf::base::size_type level_coarse)
 Interpolate a vector of DOFs on a finer mesh. More...
 
template<typename SCALAR >
std::ostream & operator<< (std::ostream &o, const ScalarFESpace< SCALAR > &fes)
 output operator for scalar parametric finite element space More...
 

Variables

static const unsigned int kout_l2_qr = 1
 
static const unsigned int kout_l2_rsfvals = 2
 
std::shared_ptr< spdlog::logger > diffusion_element_matrix_provider_logger
 logger for DiffusionElementMatrixProvider More...
 
std::shared_ptr< spdlog::logger > mass_element_matrix_provider_logger
 logger for MassElementMatrixProvider More...
 
std::shared_ptr< spdlog::logger > mass_edge_matrix_provider_logger
 logger for MassEdgeMatrixProvider More...
 
std::shared_ptr< spdlog::logger > scalar_load_element_vector_provider_logger
 logger used by ScalarLoadElementVectorProvider More...
 
std::shared_ptr< spdlog::logger > scalar_load_edge_vector_provider_logger
 logger for ScalarLoadEdgeVectorProvider class template. More...
 

Detailed Description

Collects data structures and algorithms designed for scalar finite element methods primarily meant for second-order elliptic boundary value problems.

This namespace contains a number of classes/functions which can be used to solve boundary vlaue problems with scalar finite elements. This means that the shape functions must be scalar valued, but the shape functions of a given approximation space may depend on the location in the mesh instead of only the corresponding reference element. The lf::uscalfe namespace is a specialization of this namespace to uniform scalar finite elements. (Mainly Lagrangian FE).

Examples of approximation spaces that the methodsclasses in this namespace can represent/handle are:

Typedef Documentation

◆ dim_t

Type for (co-)dimensions

Definition at line 56 of file fe.h.

◆ gdof_idx_t

Type for indices into global matrices/vectors

Definition at line 50 of file fe.h.

◆ glb_idx_t

Type for global index of entities

Definition at line 58 of file fe.h.

◆ ldof_idx_t

Type for indices referring to entity matrices/vectors

Definition at line 52 of file fe.h.

◆ size_type

Type for vector length/matrix sizes

Definition at line 54 of file fe.h.

◆ sub_idx_t

Type for indexing sub-entities

Definition at line 60 of file fe.h.

Function Documentation

◆ DiffusionElementMatrixProvider()

template<class PTR , class DIFF_COEFF >
lf::fe::DiffusionElementMatrixProvider ( PTR  fe_space,
DIFF_COEFF  alpha 
) -> DiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF >

◆ ijacobi()

double lf::fe::ijacobi ( unsigned  n,
double  alpha,
double  x 
)

Evaluate the integral of the (n-1)-th degree Jacobi Polynomial for \( *\beta = 0 \).

Parameters
nThe degree of the integrated polynomial
alphaThe \( \alpha \) parameter of the Jacobi polynomial
xThe evaluation coordinate

The integral is evaluated using

\[ \begin{aligned} L_1^\alpha(x) &= x \\ L_p^\alpha(x) &= a_pP_p^\alpha(x) + b_pP_{p-1}^\alpha(x) - *c_pP_{p-2}^\alpha(x) \end{aligned} \]

where the coefficients are defined as

\[ \begin{aligned} a_p &= \frac{p+\alpha}{(2p+\alpha-1)(2p+\alpha)} \\ b_p &= \frac{\alpha}{(2p+\alpha-1)(2p+\alpha)} \\ c_p &= \frac{p-1}{(2p+\alpha-2)(2p+\alpha-1)} \end{aligned} \]

Definition at line 206 of file hierarchic_fe.cc.

Referenced by lf::fe::FeHierarchicTria< SCALAR >::EvalReferenceShapeFunctions(), and lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions().

◆ ijacobi_dx()

double lf::fe::ijacobi_dx ( unsigned  n,
double  alpha,
double  x 
)

Computes the derivative of the n-th integrated scaled Jacobi polynomial.

Parameters
nThe degree of the integrated scaled Jacobi polynomial
alphaThe \( \alpha \) parameter of the Jacobi polynomial
xThe evaluation coordinate

The derivative is simply given by \( \frac{\partial}{\partial x} L_n^{(\alpha,0)}(x) = P_{n-1}^{(\alpha,0)}(x) \)

Definition at line 253 of file hierarchic_fe.cc.

References jacobi().

Referenced by lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions(), and lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().

◆ ilegendre()

double lf::fe::ilegendre ( unsigned  n,
double  x,
double  t 
)

computes the integral of the (n-1)-th degree scaled Legendre Polynomial

Parameters
nThe degree of the integrated polynomial
xThe evaluation coordinate
tThe scaling parameter

The integral is evaluated using

\[ \begin{aligned} L_1(x) &= x \\ 2(2n-1)L_n(x) &= P_n(x) - t^2P_{n-2}(x) \end{aligned} \]

Definition at line 60 of file hierarchic_fe.cc.

References legendre().

Referenced by lf::fe::FeHierarchicSegment< SCALAR >::EvalReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::EvalReferenceShapeFunctions(), and lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions().

◆ ilegendre_dt()

double lf::fe::ilegendre_dt ( unsigned  n,
double  x,
double  t 
)

Computes \( \frac{\partial}{\partial t} L(x;t) \).

Parameters
nThe degree of the integrated scaled Legendre polynomial
xThe evaluation coordinate
tThe scaling parameter

The derivative is given by

\[ \begin{aligned} \frac{\partial}{\partial t} L_1(x;t) &= 0 \\ \frac{\partial}{\partial t} L_n(x;t) &= -\frac{1}{2} \left( *P_{n-1}(x;t) + tP_{n-2}(x;t) \right) \end{aligned} \]

Definition at line 94 of file hierarchic_fe.cc.

References legendre().

Referenced by lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions().

◆ ilegendre_dx()

double lf::fe::ilegendre_dx ( unsigned  n,
double  x,
double  t 
)

Computes \( \frac{\partial}{\partial x} L(x;t) \).

Parameters
nThe degree of the integrated scaled Legendre polynomial
xThe evaluation coordinate
tThe scaling parameter

The derivative is simply given by \( \frac{\partial}{\partial x} L_n(x;t) = P_{n-1}(x;t) \)

Definition at line 77 of file hierarchic_fe.cc.

References legendre().

Referenced by lf::fe::FeHierarchicSegment< SCALAR >::GradientsReferenceShapeFunctions(), lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions(), lf::fe::FeHierarchicSegment< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicQuad< SCALAR >::NodalValuesToDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs().

◆ InitEssentialConditionFromFunction()

template<typename SCALAR , typename EDGESELECTOR , typename FUNCTION >
std::vector< std::pair< bool, SCALAR > > lf::fe::InitEssentialConditionFromFunction ( const lf::fe::ScalarFESpace< SCALAR > &  fes,
EDGESELECTOR &&  esscondflag,
FUNCTION &&  g 
)

Initialization of flags/values for dofs of a Scalar finite element space whose values are imposed by a specified function.

Template Parameters
SCALARscalar type of the basis functions of the finite element space
EDGESELECTORpredicate returning true for edges with fixed dofs
FUNCTIONMeshFunction which defines the imposed values on the edges
Parameters
fesThe ScalarFESpace whose dofs should be fixed.
esscondflagpredicate object whose evaluation operator returns true for all edges whose associated degrees of freedom should be set a fixed value.
gA scalar valued MeshFunction which describes the values on the edges to which the dofs should be fixed.
Returns
a vector of flag-value pairs, a true first component indicating a fixed dof, with the second component providing the value in this case.

This function interpolates a scalar-valued function into a Scalar finite element space restricted to a set of edges. It relies on the method ScalarReferenceFiniteElement::NodalValuesToDofs().

The main use of this function is the interpolation of Dirichlet data on the Dirichlet part of the boundary of a domain.

Template parameter type requirements

  • SCALAR must be a type like complex<double>
  • EDGESELECTOR must be compatible with std::function<bool(const Entity &)>
  • FUNCTION is a scalar valued MeshFunction which is evaluated on edges

This function is meant to supply the information needed for the elimination of Dirichlet boundary conditions by means of the function lf::assemble::FixFlaggedSolutionComponents().

Example

auto mesh_factory = std::make_unique<mesh::hybrid2d::MeshFactory>(2);
auto gmsh_reader = io::GmshReader(std::move(mesh_factory), "mesh.msh");
auto mesh = gmsh_reader.mesh();
auto fe_space =
std::make_shared<lf::uscalfe::FeSpaceLagrangeO1<double>>(mesh);
// We want to solve the pde
// - \laplace u = 0
// u = cos(x)*sin(y) on the boundary
// 1) Setup a mesh function that represents the prescribed values of u on the
// boundary
auto mf_boundary =
mesh::utils::MeshFunctionGlobal([&](const Eigen::Vector2d& x) {
return std::cos(x[0]) * std::sin(x[1]);
});
// 2) determine the dofs on the boundary and to what value the should be set
auto boundary_cond = InitEssentialConditionFromFunction(
*fe_space, base::PredicateTrue{}, mf_boundary);
// 3) Assemble the stiffness matrix:
assemble::COOMatrix<double> lhs(fe_space->LocGlobMap().NumDofs(),
fe_space->LocGlobMap().NumDofs());
auto mf_one = mesh::utils::MeshFunctionConstant(1.);
auto mf_zero = mesh::utils::MeshFunctionConstant(0.);
fe_space, mf_one, mf_zero);
assemble::AssembleMatrixLocally(0, fe_space->LocGlobMap(),
fe_space->LocGlobMap(), matrix_provider, lhs);
// 4) Modify the system of equations such that the boundary values are
// enforced
Eigen::VectorXd rhs(fe_space->LocGlobMap().NumDofs());
[&](unsigned int idx) { return boundary_cond[idx]; }, lhs, rhs);
// 5) solve the problem:
auto lhs_sparse = lhs.makeSparse();
Eigen::SimplicialLDLT<decltype(lhs_sparse)> solver;
solver.compute(lhs_sparse);
auto x = solver.solve(rhs);
A Function Object that can be invoked with any arguments and that always returns the value true.
Reads a Gmsh *.msh file into a mesh::MeshFactory and provides a link between mesh::Entity objects and...
Definition: gmsh_reader.h:70
MeshFunction wrapper for a simple function of physical coordinates.
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 FixFlaggedSolutionComponents(SELECTOR &&selectvals, COOMatrix< SCALAR > &A, RHSVECTOR &b)
enforce prescribed solution components
Definition: fix_dof.h:87
std::vector< std::pair< bool, SCALAR > > InitEssentialConditionFromFunction(const lf::fe::ScalarFESpace< SCALAR > &fes, EDGESELECTOR &&esscondflag, FUNCTION &&g)
Initialization of flags/values for dofs of a Scalar finite element space whose values are imposed by ...
Definition: fe_tools.h:301
ReactionDiffusionElementMatrixProvider(PTR fe_space, DIFF_COEFF alpha, REACTION_COEFF gamma) -> ReactionDiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF, REACTION_COEFF >

Definition at line 301 of file fe_tools.h.

References lf::mesh::Mesh::Entities(), lf::assemble::DofHandler::GlobalDofIndices(), lf::fe::ScalarFESpace< SCALAR >::LocGlobMap(), lf::fe::ScalarFESpace< SCALAR >::Mesh(), lf::assemble::DofHandler::NumDofs(), lf::assemble::DofHandler::NumLocalDofs(), and lf::fe::ScalarFESpace< SCALAR >::ShapeFunctionLayout().

Referenced by projects::dpg::InitEssentialConditionsFromFunctions().

◆ IntegrateMeshFunction() [1/2]

template<class MF , class QR_SELECTOR , class ENTITY_PREDICATE = base::PredicateTrue, class = std::enable_if_t< std::is_invocable_v<QR_SELECTOR, const mesh::Entity &>>>
auto lf::fe::IntegrateMeshFunction ( const lf::mesh::Mesh mesh,
const MF &  mf,
const QR_SELECTOR &  qr_selector,
const ENTITY_PREDICATE &  ep = base::PredicateTrue{},
int  codim = 0 
) -> mesh::utils::MeshFunctionReturnType<MF>

Integrate a MeshFunction over a mesh (with quadrature rules)

Template Parameters
MFThe type of the mesh function.
QR_SELECTORThe type of qr_selector (see below)
ENTITY_PREDICATEThe type of the entity predicate (see below)
Parameters
meshThe mesh to integrate over
mfThe mesh function to integrate
qr_selectorProvides the quadrature rule for every entity of the mesh.
epSelects the entities over which mf is integrated (default: all entities)
codimThe codimension of the entities over which mf is integrated.
Returns
The integrated value

Requirements for QR_SELECTOR

QR_SELECTOR should overload operator() as follows:

quad::QuadRule operator()(const mesh::Entity& e) const
Interface class representing a topological entity in a cellular complex
Definition: entity.h:39
Represents a Quadrature Rule over one of the Reference Elements.
Definition: quad_rule.h:58

i.e. it should return the quadrature rule for every entity e of the mesh that is to be used for computing the integral of mf over e.

Requirements for ENTITY_PREDICATE

The entity predicate should overload operator() as follows:

bool operator()(const mesh::Entity& e) const

It should return true, if e is part of the integration domain and false if it is not.

Example

auto mesh_factory = std::make_unique<mesh::hybrid2d::MeshFactory>(2);
auto gmsh_reader = io::GmshReader(std::move(mesh_factory), "mesh.msh");
auto mesh = gmsh_reader.mesh();
// integrate the function sin(x)*cos(y) over the mesh using 5th-degree
// quadrature rules
[](const Eigen::Vector2d& x) { return std::sin(x[0]) * std::cos(x[1]); });
// select the quadrature rule explicitly for every element:
auto integral = IntegrateMeshFunction(*mesh, mf, [](const mesh::Entity& e) {
return quad::make_QuadRule(e.RefEl(), 5);
});
virtual base::RefEl RefEl() const =0
Describes the reference element type of this entity.
auto IntegrateMeshFunction(const lf::mesh::Mesh &mesh, const MF &mf, const QR_SELECTOR &qr_selector, const ENTITY_PREDICATE &ep=base::PredicateTrue{}, int codim=0) -> mesh::utils::MeshFunctionReturnType< MF >
Integrate a MeshFunction over a mesh (with quadrature rules)
Definition: fe_tools.h:111
QuadRule make_QuadRule(base::RefEl ref_el, unsigned degree)
Returns a QuadRule object for the given Reference Element and Degree.

Definition at line 111 of file fe_tools.h.

Referenced by projects::ipdg_stokes::post_processing::DGnorm(), projects::ipdg_stokes::post_processing::L2norm(), projects::dpg::test::TestConververgencePrimalDPGAdaptedNormConvectionDiffusionDirichletBVP(), and projects::dpg::test::TestConververgencePrimalDPGConvectionDiffusionDirichletBVP().

◆ IntegrateMeshFunction() [2/2]

template<class MF , class ENTITY_PREDICATE = base::PredicateTrue>
auto lf::fe::IntegrateMeshFunction ( const lf::mesh::Mesh mesh,
const MF &  mf,
int  quad_degree,
const ENTITY_PREDICATE &  ep = base::PredicateTrue{},
int  codim = 0 
)

Integrate a mesh function over a mesh using quadrature rules of uniform order.

Template Parameters
MFtype of mesh function to integrate
ENTITY_PREDICATEtype of entity predicate (see below)
Parameters
meshThe mesh over which mf is integrated.
mfThe mesh function which is integrated
quad_degreeThe quadrature degree of the quadrature rules that are to be used for integration. Internally Gauss-rules created by quad::make_QuadRule are used.
epThe entity predicate selecting the entities over which mf is integrated.
codimThe codimension of the entities over which mf is integrated.
Returns
mf integrated over the entities e of mf where ep(e)==true.

Requirements for ENTITY_PREDICATE

The entity predicate should overload operator() as follows:

bool operator()(const mesh::Entity& e) const

It should return true, if e is part of the integration domain and false if it is not.

Example

auto mesh_factory = std::make_unique<mesh::hybrid2d::MeshFactory>(2);
auto gmsh_reader = io::GmshReader(std::move(mesh_factory), "mesh.msh");
auto mesh = gmsh_reader.mesh();
// integrate the function sin(x)*cos(y) over the mesh using 5th-degree
// quadrature rules
[](const Eigen::Vector2d& x) { return std::sin(x[0]) * std::cos(x[1]); });
auto integral = IntegrateMeshFunction(*mesh, mf, 5);

Definition at line 157 of file fe_tools.h.

◆ jacobi() [1/2]

double lf::fe::jacobi ( unsigned  n,
double  alpha,
double  beta,
double  x 
)

Computes the n-th degree shifted Jacobi polynomial.

Parameters
nThe degree of the polynomial
alphaThe \( \alpha \) parameter of the Jacobi polynomial
betaThe \( \beta \) parameter of the Jacobi polynomial
xThe evaluation coordinate

We use the recurrence relation for non-shifted Jacobi polynomials

\[ \begin{aligned} P_0^{(\alpha,\beta)}(x) &= 1 \\ P_1^{(\alpha,\beta)}(x) &= \frac{1}{2} \left( \alpha - \beta + *(\alpha + \beta + 2)x \right) \\ P_{n+1}^{(\alpha,\beta)}(x) &= \frac{1}{a_n} *\left( (b_n+c_nx)P_n^{(\alpha,\beta)}(x) - d_nP_{n-1}^{(\alpha,\beta)}(x) *\right) \end{aligned} \]

where

\[ \begin{aligned} a_n &= 2(n+1)(n+\alpha+\beta+1)(2n+\alpha+\beta) \\ b_n &= (2n+\alpha+\beta+1)(\alpha^2-\beta^2) \\ c_n &= (2n+\alpha+\beta)(2n+\alpha+\beta+1)(2n+\alpha+\beta+2) \\ d_n &= 2(n+\alpha)(n+\beta)(2n+\alpha+\beta+2) \end{aligned} \]

Definition at line 145 of file hierarchic_fe.cc.

Referenced by ijacobi_dx(), jacobi(), and jacobi_dx().

◆ jacobi() [2/2]

double lf::fe::jacobi ( unsigned  n,
double  alpha,
double  x 
)

Computes the n-th degree shifted Jacobi polynomial for \( \beta = 0 \).

Parameters
nThe degree of the polynomial
alphaThe \( \alpha \) parameter of the Jacobi polynomial
xThe evaluation coordinate

Definition at line 181 of file hierarchic_fe.cc.

References jacobi().

◆ jacobi_dx()

double lf::fe::jacobi_dx ( unsigned  n,
double  alpha,
double  x 
)

Computes the derivative of the n-th degree Jacobi Polynomial for \( *\beta = 0 \).

Parameters
nThe degree of the differentiated polynomial
alphaThe \( \alpha \) parameter of the Jacobi Polynomial
xThe evaluation coordinate

The derivative is evaluated using

\[ {P^{(\alpha,0)}_n}'(x) = \frac{\alpha+n+1}{2} P^{(\alpha+1,1)}_{n-1}(x) \]

Definition at line 269 of file hierarchic_fe.cc.

References jacobi().

Referenced by lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().

◆ legendre()

double lf::fe::legendre ( unsigned  n,
double  x,
double  t 
)

computes the n-th degree scaled Legendre Polynomial \( P_n(x;t) *\)

Parameters
nThe degree of the polynomial
xThe evaluation coordinate
tThe scaling parameter

To evaluate the scaled Legendre Polynomials \( P_n(x;t) \), we use that

\[ \begin{aligned} P_0(x;t) &= 1 \\ P_1(x;t) &= 2x - t \\ nP_n(x;t) &= (2n-1)(2x-t)P_{n-1}(x;t) - (n-1)t^2P_{p-2}(x;t) \end{aligned} \]

Definition at line 27 of file hierarchic_fe.cc.

Referenced by ilegendre(), ilegendre_dt(), ilegendre_dx(), and legendre_dx().

◆ legendre_dx()

double lf::fe::legendre_dx ( unsigned  n,
double  x,
double  t 
)

Computes the derivative of the n-th degree scaled Legendre polynomial.

Parameters
nThe degree of the polynomial
xThe evaluation coordinate
tThe scaling parameter

The derivative is given by

\[ \begin{aligned} \frac{\partial}{\partial x} L_0(x;t) &= 0 \\ \frac{\partial}{\partial x} L_n(x;t) &= 2nP_{n-1}(x;t) + *(2x-t)\frac{\partial}{\partial x}P_{n-1}(x;t) \\ \end{aligned} \]

Definition at line 114 of file hierarchic_fe.cc.

References legendre(), and legendre_dx().

Referenced by legendre_dx(), lf::fe::FeHierarchicSegment< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicQuad< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().

◆ MassEdgeMatrixProvider()

template<class PTR , class COEFF , class EDGESELECTOR = base::PredicateTrue>
lf::fe::MassEdgeMatrixProvider ( PTR  ,
COEFF  coeff,
EDGESELECTOR  edge_predicate = base::PredicateTrue{} 
) -> MassEdgeMatrixProvider< typename PTR::element_type::Scalar, COEFF, EDGESELECTOR >

◆ MassElementMatrixProvider()

template<class PTR , class REACTION_COEFF >
lf::fe::MassElementMatrixProvider ( PTR  fe_space,
REACTION_COEFF  gamma 
) -> MassElementMatrixProvider< typename PTR::element_type::Scalar, REACTION_COEFF >

◆ MeshFunctionFE()

template<class T , class SCALAR_COEFF >
lf::fe::MeshFunctionFE ( std::shared_ptr< T >  ,
const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &   
) -> MeshFunctionFE< typename T::Scalar, SCALAR_COEFF >

◆ MeshFunctionGradFE()

template<class T , class SCALAR_COEFF >
lf::fe::MeshFunctionGradFE ( std::shared_ptr< T >  ,
const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &   
) -> MeshFunctionGradFE< typename T::Scalar, SCALAR_COEFF >

◆ NodalProjection()

template<typename SCALAR , typename MF , typename SELECTOR = base::PredicateTrue>
auto lf::fe::NodalProjection ( const lf::fe::ScalarFESpace< SCALAR > &  fe_space,
MF &&  u,
SELECTOR &&  pred = base::PredicateTrue{} 
)

Computes nodal projection of a mesh function and returns the finite element basis expansion coefficients of the result.

Template Parameters
SCALARa scalar type
MFa MeshFunction representing the scalar valued function that should be projected
SELECTORpredicate type for selecting cells to be visited
Parameters
fe_spacea uniform Lagrangian finite element space, providing finite element specifications for the cells of the mesh
ufunctor object supplying a scalar-valued function that is to be projected
predpredicate object for the selection of relevant cells
Returns
column vector of basis expansion coefficients for the resulting finite element function

The implementation relies on the method ScalarReferenceFiniteElement::NodalValuesToDofs(). Refer to its documentation. This method is called for each active cell to set the coefficients for the global shape functions associated with that cell.

Example

auto mesh_factory = std::make_unique<mesh::hybrid2d::MeshFactory>(2);
auto gmsh_reader = io::GmshReader(std::move(mesh_factory), "mesh.msh");
auto mesh = gmsh_reader.mesh();
// project a first-degree polynomial onto a first order lagrange space and
// make sure the representation is exact:
auto fe_space =
std::make_shared<lf::uscalfe::FeSpaceLagrangeO1<double>>(mesh);
[](const Eigen::Vector2d& x) { return 2 + 3 * x[0] + 4 * x[1]; });
auto dof_vector = NodalProjection(*fe_space, mf_linear);
auto mf_fe = lf::fe::MeshFunctionFE(fe_space, dof_vector);
assert(IntegrateMeshFunction(*mesh, squaredNorm(mf_fe - mf_linear), 2) <
1e-12);
auto NodalProjection(const lf::fe::ScalarFESpace< SCALAR > &fe_space, MF &&u, SELECTOR &&pred=base::PredicateTrue{})
Computes nodal projection of a mesh function and returns the finite element basis expansion coefficie...
Definition: fe_tools.h:199
MeshFunctionFE(std::shared_ptr< T >, const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &) -> MeshFunctionFE< typename T::Scalar, SCALAR_COEFF >

Definition at line 199 of file fe_tools.h.

Referenced by prolongate().

◆ operator<<()

template<typename SCALAR >
std::ostream & lf::fe::operator<< ( std::ostream &  o,
const ScalarFESpace< SCALAR > &  fes 
)

output operator for scalar parametric finite element space

Definition at line 100 of file scalar_fe_space.h.

◆ prolongate()

template<typename SCALAR_COEFF , typename FES_COARSE , typename FES_FINE >
Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > lf::fe::prolongate ( const lf::refinement::MeshHierarchy mh,
std::shared_ptr< FES_COARSE >  fespace_coarse,
std::shared_ptr< FES_FINE >  fespace_fine,
const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &  dofs_coarse,
lf::base::size_type  level_coarse 
)

Interpolate a vector of DOFs on a finer mesh.

Template Parameters
SCALAR_COEFFThe scalar of the coefficient vector
FES_COARSEThe FE space on the coarse mesh
FES_FINEThe FE space on the fine mesh
Parameters
mhA reference to the MeshHierarchy containing the underlying meshes
fespace_coarseThe FE space on the coarse mesh
fespace_fineThe FE space on the fine mesh
dofs_coarseA basis function coefficient vector on the coarse mesh
level_coarseThe level of the coarse mesh
Returns
An interpolated vector of basis function coefficients on the fine mesh

Objects of type refinement::MeshHierarchy contain sequences of nested meshes. A finite-element function on a coarse mesh, which can be regarded as just another continuous function on the meshed domain, can be interpolated to yield a finite element function on the next finer mesh. This function realizes the conversion of the corresponding basis expansion coefficient vectors.

Definition at line 33 of file prolongation.h.

References NodalProjection(), lf::assemble::DofHandler::NumDofs(), and lf::refinement::MeshHierarchy::NumLevels().

◆ ScalarLoadEdgeVectorProvider()

template<class PTR , class FUNCTOR , class EDGESELECTOR = base::PredicateTrue>
lf::fe::ScalarLoadEdgeVectorProvider ( PTR  ,
FUNCTOR  ,
EDGESELECTOR  = base::PredicateTrue{} 
) -> ScalarLoadEdgeVectorProvider< typename PTR::element_type::Scalar, FUNCTOR, EDGESELECTOR >

◆ ScalarLoadElementVectorProvider()

template<class PTR , class MESH_FUNCTION >
lf::fe::ScalarLoadElementVectorProvider ( PTR  fe_space,
MESH_FUNCTION  mf 
) -> ScalarLoadElementVectorProvider< typename PTR::element_type::Scalar, MESH_FUNCTION >

Variable Documentation

◆ diffusion_element_matrix_provider_logger

std::shared_ptr< spdlog::logger > lf::fe::diffusion_element_matrix_provider_logger
Initial value:
=
base::InitLogger("lf::fe::diffusion_element_matrix_provider_logger")

logger for DiffusionElementMatrixProvider

Definition at line 20 of file loc_comp_ellbvp.cc.

Referenced by lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval().

◆ kout_l2_qr

const unsigned int lf::fe::kout_l2_qr = 1
static

Definition at line 26 of file fe_tools.h.

◆ kout_l2_rsfvals

const unsigned int lf::fe::kout_l2_rsfvals = 2
static

Definition at line 27 of file fe_tools.h.

◆ mass_edge_matrix_provider_logger

std::shared_ptr< spdlog::logger > lf::fe::mass_edge_matrix_provider_logger
Initial value:
=
base::InitLogger("lf::fe::mass_edge_matrix_provider_logger")

logger for MassEdgeMatrixProvider

Definition at line 26 of file loc_comp_ellbvp.cc.

◆ mass_element_matrix_provider_logger

std::shared_ptr< spdlog::logger > lf::fe::mass_element_matrix_provider_logger
Initial value:
=
base::InitLogger("lf::fe::mass_element_matrix_provider_logger")

logger for MassElementMatrixProvider

Definition at line 23 of file loc_comp_ellbvp.cc.

Referenced by lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval().

◆ scalar_load_edge_vector_provider_logger

std::shared_ptr< spdlog::logger > lf::fe::scalar_load_edge_vector_provider_logger
Initial value:
=
base::InitLogger("lf::fe::scalar_load_edge_vector_provider_logger")

logger for ScalarLoadEdgeVectorProvider class template.

Definition at line 32 of file loc_comp_ellbvp.cc.

◆ scalar_load_element_vector_provider_logger

std::shared_ptr< spdlog::logger > lf::fe::scalar_load_element_vector_provider_logger
Initial value:
=
base::InitLogger("lf::fe::scalar_load_element_vector_provider_logger")

logger used by ScalarLoadElementVectorProvider

Definition at line 29 of file loc_comp_ellbvp.cc.

Referenced by lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval().