LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
Public Types | Public Member Functions | Private Attributes | List of all members
lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF > Class Template Reference

A MeshFunction representing an element from a ScalarFESpace (e.g. solution of BVP) More...

#include <lf/fe/fe.h>

Public Types

using Scalar = decltype(SCALAR_FE(0) * SCALAR_COEFF(0))
 

Public Member Functions

 MeshFunctionFE (std::shared_ptr< const ScalarFESpace< SCALAR_FE > > fe_space, Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > coeff_vector)
 Create a new MeshFunctionFE from a ScalarFESpace and a coefficient vector.
 
std::vector< Scalaroperator() (const lf::mesh::Entity &e, const Eigen::MatrixXd &local) const
 Evaluate the mesh function on a MeshEntity.
 
std::shared_ptr< const mesh::MeshgetMesh () const
 Convenience method to retrieve the underlying mesh.
 
std::shared_ptr< const ScalarFESpace< SCALAR_FE > > getFESpace () const
 Convenience method to retrieve the finite element space in which the mesh function lives.
 

Private Attributes

std::shared_ptr< const ScalarFESpace< SCALAR_FE > > fe_space_
 
Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > dof_vector_
 

Detailed Description

template<class SCALAR_FE, class SCALAR_COEFF>
class lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >

A MeshFunction representing an element from a ScalarFESpace (e.g. solution of BVP)

Template Parameters
SCALAR_FEThe scalar type of the finite element basis functions.
SCALAR_COEFFThe scalar type of the coefficient vector

The MeshFunctionFE takes essentially two parameters:

Note
A MeshFunctionFE can be evaluated on entities of all codimensions. The mesh function will simply use the ScalarReferenceFiniteElement from the feSpace.

Use case

double computeL2ErrorNorm(
std::shared_ptr<const lf::uscalfe::UniformScalarFESpace<double>> fe_space_p,
const Eigen::VectorXd& mu) {
// Reference to underlying mesh
const lf::mesh::Mesh& mesh{*(fe_space_p->Mesh())};
// Lambda function providing the "exact solution"
auto u = [](Eigen::Vector2d x) -> double {
return std::log(x[0] * x[0] + x[1] + 1.0);
};
// Has to be wrapped into a mesh function for error computation
// create mesh functions representing solution
auto mf_sol = lf::fe::MeshFunctionFE<double, double>(fe_space_p, mu);
// compute errors with 10-th order quadrature rules
double L2err_2 = // NOLINT
std::sqrt(IntegrateMeshFunction(mesh, squaredNorm(mf_sol - mf_u), 2));
return std::sqrt(L2err_2);
}

Definition at line 44 of file mesh_function_fe.h.

Member Typedef Documentation

◆ Scalar

template<class SCALAR_FE , class SCALAR_COEFF >
using lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::Scalar = decltype(SCALAR_FE(0) * SCALAR_COEFF(0))

Definition at line 48 of file mesh_function_fe.h.

Constructor & Destructor Documentation

◆ MeshFunctionFE()

template<class SCALAR_FE , class SCALAR_COEFF >
lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::MeshFunctionFE ( std::shared_ptr< const ScalarFESpace< SCALAR_FE > > fe_space,
Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > coeff_vector )
inline

Create a new MeshFunctionFE from a ScalarFESpace and a coefficient vector.

Parameters
fe_spacethe approximation space in which the function lies.
coeff_vectorDefines the coefficients in front of the basis functions of fe_space
Warning
This constructor just takes a reference to the vector of basis expansion coefficients. Thus, this vector has to be "kept alive" as long as the mesh function exists.

Definition at line 61 of file mesh_function_fe.h.

Member Function Documentation

◆ getFESpace()

template<class SCALAR_FE , class SCALAR_COEFF >
std::shared_ptr< const ScalarFESpace< SCALAR_FE > > lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::getFESpace ( ) const
inline

Convenience method to retrieve the finite element space in which the mesh function lives.

Returns
shared_ptr to ScalarFESpace in which the mesh function lives.

Definition at line 107 of file mesh_function_fe.h.

References lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::fe_space_.

◆ getMesh()

template<class SCALAR_FE , class SCALAR_COEFF >
std::shared_ptr< const mesh::Mesh > lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::getMesh ( ) const
inline

Convenience method to retrieve the underlying mesh.

Returns
The mesh on which this mesh function is defined.

Definition at line 97 of file mesh_function_fe.h.

References lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::fe_space_.

◆ operator()()

template<class SCALAR_FE , class SCALAR_COEFF >
std::vector< Scalar > lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::operator() ( const lf::mesh::Entity & e,
const Eigen::MatrixXd & local ) const
inline

Evaluate the mesh function on a MeshEntity.

Parameters
ethe relevant mesh entity
localreference coordinates of evalation points passed in the columns of a matrix
Returns
values of the function at evaluation points

Definition at line 72 of file mesh_function_fe.h.

References lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::dof_vector_, and lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::fe_space_.

Member Data Documentation

◆ dof_vector_

template<class SCALAR_FE , class SCALAR_COEFF >
Eigen::Matrix<SCALAR_COEFF, Eigen::Dynamic, 1> lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::dof_vector_
private

◆ fe_space_

template<class SCALAR_FE , class SCALAR_COEFF >
std::shared_ptr<const ScalarFESpace<SCALAR_FE> > lf::fe::MeshFunctionFE< SCALAR_FE, SCALAR_COEFF >::fe_space_
private

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