LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
|
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 |
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 ScalarFESpace (e.g. solution of BVP) More... | |
class | MeshFunctionGradFE |
A MeshFunction representing the gradient of a function from a scalar finite element space (e.g. gradient of a solution of BVP). More... | |
class | ScalarFESpace |
Space of scalar valued finite element functions on a 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<mesh::utils::MeshFunction MF, class QR_SELECTOR , class ENTITY_PREDICATE = base::PredicateTrue> requires 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) | |
template<mesh::utils::MeshFunction 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) -> mesh::utils::MeshFunctionReturnType< MF > |
Integrate a mesh function over a mesh using quadrature rules of uniform order. | |
template<typename SCALAR , mesh::utils::MeshFunction MF, typename SELECTOR = base::PredicateTrue> | |
auto | NodalProjection (const lf::fe::ScalarFESpace< SCALAR > &fe_space, MF const &u, SELECTOR &&pred=base::PredicateTrue{}) -> Eigen::Vector< decltype(SCALAR{0} *mesh::utils::MeshFunctionReturnType< std::remove_reference_t< MF > >{0}), Eigen::Dynamic > |
Computes nodal projection of a mesh function and returns the finite element basis expansion coefficients of the result. | |
template<typename SCALAR , typename EDGESELECTOR , mesh::utils::MeshFunction FUNCTION> | |
std::vector< std::pair< bool, SCALAR > > | InitEssentialConditionFromFunction (const lf::fe::ScalarFESpace< SCALAR > &fes, EDGESELECTOR &&esscondflag, FUNCTION const &g) |
Initialization of flags/values for dofs of a Scalar finite element space whose values are imposed by a specified function. | |
double | Legendre (unsigned n, double x, double t=1) |
computes the n -th degree scaled Legendre Polynomial \( P_n(x;t)
\) | |
double | ILegendre (unsigned n, double x, double t=1) |
computes the integral of the (n-1)-th degree scaled Legendre Polynomial | |
double | ILegendreDx (unsigned n, double x, double t=1) |
Computes \( \frac{\partial}{\partial x} L_n(x;t) \). | |
double | ILegendreDt (unsigned n, double x, double t=1) |
Computes \( \frac{\partial}{\partial t} L_n(x;t) \). | |
double | LegendreDx (unsigned n, double x, double t=1) |
Computes the derivative of the n-th degree scaled Legendre polynomial. | |
double | Jacobi (unsigned n, double alpha, double beta, double x) |
Computes the n-th degree shifted Jacobi polynomial. | |
double | Jacobi (unsigned n, double alpha, double x) |
Computes the n-th degree shifted Jacobi polynomial for \( \beta = 0
\). | |
double | IJacobi (unsigned n, double alpha, double x) |
Evaluate the integral of the (n-1)-th degree Jacobi Polynomial for \(
\beta = 0 \). | |
double | IJacobiDx (unsigned n, double alpha, double x) |
Computes the derivative of the n-th integrated scaled Jacobi polynomial. | |
double | JacobiDx (unsigned n, double alpha, double x) |
Computes the derivative of the n-th degree Jacobi Polynomial for \(
\beta = 0 \). | |
std::shared_ptr< spdlog::logger > & | DiffusionElementMatrixProviderLogger () |
logger for DiffusionElementMatrixProvider | |
std::shared_ptr< spdlog::logger > & | MassElementMatrixProviderLogger () |
logger for MassElementMatrixProvider | |
std::shared_ptr< spdlog::logger > & | MassEdgeMatrixProviderLogger () |
logger for MassEdgeMatrixProvider | |
std::shared_ptr< spdlog::logger > & | ScalarLoadElementVectorProviderLogger () |
logger used by ScalarLoadElementVectorProvider | |
std::shared_ptr< spdlog::logger > & | ScalarLoadEdgeVectorProviderLogger () |
logger for ScalarLoadEdgeVectorProvider class template. | |
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. | |
template<typename SCALAR > | |
std::ostream & | operator<< (std::ostream &o, const ScalarFESpace< SCALAR > &fes) |
output operator for scalar parametric finite element space | |
Variables | |
static const unsigned int | kout_l2_qr = 1 |
static const unsigned int | kout_l2_rsfvals = 2 |
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:
lf::fe::DiffusionElementMatrixProvider | ( | PTR | fe_space, |
DIFF_COEFF | alpha ) -> DiffusionElementMatrixProvider< typename PTR::element_type::Scalar, DIFF_COEFF > |
std::shared_ptr< spdlog::logger > & lf::fe::DiffusionElementMatrixProviderLogger | ( | ) |
logger for DiffusionElementMatrixProvider
Definition at line 20 of file loc_comp_ellbvp.cc.
References lf::base::InitLogger().
Referenced by lf::fe::DiffusionElementMatrixProvider< SCALAR, DIFF_COEFF >::Eval().
double lf::fe::IJacobi | ( | unsigned | n, |
double | alpha, | ||
double | x ) |
Evaluate the integral of the (n-1)-th degree Jacobi Polynomial for \( \beta = 0 \).
n | The degree of the integrated polynomial |
alpha | The \( \alpha \) parameter of the Jacobi polynomial |
x | The 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 93 of file hierarchic_fe.cc.
Referenced by lf::fe::FeHierarchicTria< SCALAR >::EvalReferenceShapeFunctions(), and lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions().
double lf::fe::IJacobiDx | ( | unsigned | n, |
double | alpha, | ||
double | x ) |
Computes the derivative of the n-th integrated scaled Jacobi polynomial.
n | The degree of the integrated scaled Jacobi polynomial |
alpha | The \( \alpha \) parameter of the Jacobi polynomial |
x | The 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 130 of file hierarchic_fe.cc.
References Jacobi().
Referenced by lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions(), and lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().
double lf::fe::ILegendre | ( | unsigned | n, |
double | x, | ||
double | t = 1 ) |
computes the integral of the (n-1)-th degree scaled Legendre Polynomial
n | The degree of the integrated polynomial |
x | The evaluation coordinate |
t | The 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 31 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().
double lf::fe::ILegendreDt | ( | unsigned | n, |
double | x, | ||
double | t = 1 ) |
Computes \( \frac{\partial}{\partial t} L_n(x;t) \).
n | The degree of the integrated scaled Legendre polynomial |
x | The evaluation coordinate |
t | The 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 43 of file hierarchic_fe.cc.
References Legendre().
Referenced by lf::fe::FeHierarchicTria< SCALAR >::GradientsReferenceShapeFunctions().
double lf::fe::ILegendreDx | ( | unsigned | n, |
double | x, | ||
double | t = 1 ) |
Computes \( \frac{\partial}{\partial x} L_n(x;t) \).
n | The degree of the integrated scaled Legendre polynomial |
x | The evaluation coordinate |
t | The scaling parameter |
The derivative is simply given by \( \frac{\partial}{\partial x} L_n(x;t) = P_{n-1}(x;t) \)
Definition at line 39 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().
std::vector< std::pair< bool, SCALAR > > lf::fe::InitEssentialConditionFromFunction | ( | const lf::fe::ScalarFESpace< SCALAR > & | fes, |
EDGESELECTOR && | esscondflag, | ||
FUNCTION const & | g ) |
Initialization of flags/values for dofs of a Scalar finite element space whose values are imposed by a specified function.
SCALAR | scalar type of the basis functions of the finite element space |
EDGESELECTOR | predicate returning true for edges with fixed dofs |
FUNCTION | MeshFunction which defines the imposed values on the edges |
fes | The ScalarFESpace whose dofs should be fixed. |
esscondflag | predicate object whose evaluation operator returns true for all edges whose associated degrees of freedom should be set a fixed value. |
g | A scalar valued MeshFunction which describes the values on the edges to which the dofs should be fixed. |
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.
complex<double>
std::function<bool(const Entity &)>
This function is meant to supply the information needed for the elimination of Dirichlet boundary conditions by means of the function lf::assemble::FixFlaggedSolutionComponents().
Definition at line 303 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().
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)
MF | The type of the mesh function. |
QR_SELECTOR | The type of qr_selector (see below) |
ENTITY_PREDICATE | The type of the entity predicate (see below) |
mesh | The mesh to integrate over |
mf | The mesh function to integrate |
qr_selector | Provides the quadrature rule for every entity of the mesh. |
ep | Selects the entities over which mf is integrated (default: all entities) |
codim | The codimension of the entities over which mf is integrated. |
QR_SELECTOR
should overload operator()
as follows:
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
.
The entity predicate should overload operator()
as follows:
It should return true
, if e
is part of the integration domain and false
if it is not.
Definition at line 110 of file fe_tools.h.
auto lf::fe::IntegrateMeshFunction | ( | const lf::mesh::Mesh & | mesh, |
const MF & | mf, | ||
int | quad_degree, | ||
const ENTITY_PREDICATE & | ep = base::PredicateTrue{}, | ||
int | codim = 0 ) -> mesh::utils::MeshFunctionReturnType<MF> |
Integrate a mesh function over a mesh using quadrature rules of uniform order.
MF | type of mesh function to integrate |
ENTITY_PREDICATE | type of entity predicate (see below) |
mesh | The mesh over which mf is integrated. |
mf | The mesh function which is integrated |
quad_degree | The quadrature degree of the quadrature rules that are to be used for integration. Internally Gauss-rules created by quad::make_QuadRule are used. |
ep | The entity predicate selecting the entities over which mf is integrated. |
codim | The codimension of the entities over which mf is integrated. |
mf
integrated over the entities e
of mf
where ep(e)==true
.The entity predicate should overload operator()
as follows:
It should return true
, if e
is part of the integration domain and false
if it is not.
Definition at line 155 of file fe_tools.h.
double lf::fe::Jacobi | ( | unsigned | n, |
double | alpha, | ||
double | beta, | ||
double | x ) |
Computes the n-th degree shifted Jacobi polynomial.
n | The degree of the polynomial |
alpha | The \( \alpha \) parameter of the Jacobi polynomial |
beta | The \( \beta \) parameter of the Jacobi polynomial |
x | The 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 60 of file hierarchic_fe.cc.
Referenced by IJacobiDx(), Jacobi(), and JacobiDx().
double lf::fe::Jacobi | ( | unsigned | n, |
double | alpha, | ||
double | x ) |
Computes the n-th degree shifted Jacobi polynomial for \( \beta = 0 \).
n | The degree of the polynomial |
alpha | The \( \alpha \) parameter of the Jacobi polynomial |
x | The evaluation coordinate |
Definition at line 89 of file hierarchic_fe.cc.
References Jacobi().
double lf::fe::JacobiDx | ( | unsigned | n, |
double | alpha, | ||
double | x ) |
Computes the derivative of the n-th degree Jacobi Polynomial for \( \beta = 0 \).
n | The degree of the differentiated polynomial |
alpha | The \( \alpha \) parameter of the Jacobi Polynomial |
x | The 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 134 of file hierarchic_fe.cc.
References Jacobi().
Referenced by lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().
double lf::fe::Legendre | ( | unsigned | n, |
double | x, | ||
double | t = 1 ) |
computes the n
-th degree scaled Legendre Polynomial \( P_n(x;t)
\)
n | The degree of the polynomial |
x | The evaluation coordinate |
t | The 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 13 of file hierarchic_fe.cc.
Referenced by ILegendre(), ILegendreDt(), ILegendreDx(), and LegendreDx().
double lf::fe::LegendreDx | ( | unsigned | n, |
double | x, | ||
double | t = 1 ) |
Computes the derivative of the n-th degree scaled Legendre polynomial.
n | The degree of the polynomial |
x | The evaluation coordinate |
t | The 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 51 of file hierarchic_fe.cc.
References Legendre(), and LegendreDx().
Referenced by LegendreDx(), lf::fe::FeHierarchicSegment< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicQuad< SCALAR >::NodalValuesToDofs(), lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToEdgeDofs(), and lf::fe::FeHierarchicTria< SCALAR >::NodalValuesToFaceDofs().
lf::fe::MassEdgeMatrixProvider | ( | PTR | , |
COEFF | coeff, | ||
EDGESELECTOR | edge_predicate = base::PredicateTrue{} ) -> MassEdgeMatrixProvider< typename PTR::element_type::Scalar, COEFF, EDGESELECTOR > |
std::shared_ptr< spdlog::logger > & lf::fe::MassEdgeMatrixProviderLogger | ( | ) |
logger for MassEdgeMatrixProvider
Definition at line 32 of file loc_comp_ellbvp.cc.
References lf::base::InitLogger().
lf::fe::MassElementMatrixProvider | ( | PTR | fe_space, |
REACTION_COEFF | gamma ) -> MassElementMatrixProvider< typename PTR::element_type::Scalar, REACTION_COEFF > |
std::shared_ptr< spdlog::logger > & lf::fe::MassElementMatrixProviderLogger | ( | ) |
logger for MassElementMatrixProvider
Definition at line 26 of file loc_comp_ellbvp.cc.
References lf::base::InitLogger().
Referenced by lf::fe::MassElementMatrixProvider< SCALAR, REACTION_COEFF >::Eval().
lf::fe::MeshFunctionFE | ( | std::shared_ptr< T > | , |
const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > & | ) -> MeshFunctionFE< typename T::Scalar, SCALAR_COEFF > |
lf::fe::MeshFunctionGradFE | ( | std::shared_ptr< T > | , |
const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > & | ) -> MeshFunctionGradFE< typename T::Scalar, SCALAR_COEFF > |
auto lf::fe::NodalProjection | ( | const lf::fe::ScalarFESpace< SCALAR > & | fe_space, |
MF const & | u, | ||
SELECTOR && | pred = base::PredicateTrue{} ) -> Eigen::Vector<decltype(SCALAR{0} * mesh::utils::MeshFunctionReturnType< std::remove_reference_t<MF>>{0}), Eigen::Dynamic> |
Computes nodal projection of a mesh function and returns the finite element basis expansion coefficients of the result.
SCALAR | a scalar type |
MF | a MeshFunction representing the scalar valued function that should be projected |
SELECTOR | predicate type for selecting cells to be visited |
fe_space | a scalar finite element space, providing finite element specifications for the cells of the mesh |
u | functor object supplying a scalar-valued function that is to be projected |
pred | predicate object for the selection of relevant cells |
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.
Definition at line 198 of file fe_tools.h.
Referenced by prolongate().
std::ostream & lf::fe::operator<< | ( | std::ostream & | o, |
const ScalarFESpace< SCALAR > & | fes ) |
output operator for scalar parametric finite element space
Definition at line 120 of file scalar_fe_space.h.
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.
SCALAR_COEFF | The scalar of the coefficient vector |
FES_COARSE | The FE space on the coarse mesh |
FES_FINE | The FE space on the fine mesh |
mh | A reference to the MeshHierarchy containing the underlying meshes |
fespace_coarse | The FE space on the coarse mesh |
fespace_fine | The FE space on the fine mesh |
dofs_coarse | A basis function coefficient vector on the coarse mesh |
level_coarse | The level of the coarse 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().
lf::fe::ScalarLoadEdgeVectorProvider | ( | PTR | , |
FUNCTOR | , | ||
EDGESELECTOR | = base::PredicateTrue{} ) -> ScalarLoadEdgeVectorProvider< typename PTR::element_type::Scalar, FUNCTOR, EDGESELECTOR > |
std::shared_ptr< spdlog::logger > & lf::fe::ScalarLoadEdgeVectorProviderLogger | ( | ) |
logger for ScalarLoadEdgeVectorProvider class template.
Definition at line 43 of file loc_comp_ellbvp.cc.
References lf::base::InitLogger().
lf::fe::ScalarLoadElementVectorProvider | ( | PTR | fe_space, |
MESH_FUNCTION | mf ) -> ScalarLoadElementVectorProvider< typename PTR::element_type::Scalar, MESH_FUNCTION > |
std::shared_ptr< spdlog::logger > & lf::fe::ScalarLoadElementVectorProviderLogger | ( | ) |
logger used by ScalarLoadElementVectorProvider
Definition at line 37 of file loc_comp_ellbvp.cc.
References lf::base::InitLogger().
Referenced by lf::fe::ScalarLoadElementVectorProvider< SCALAR, MESH_FUNCTION >::Eval().
|
static |
Definition at line 26 of file fe_tools.h.
|
static |
Definition at line 27 of file fe_tools.h.