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::MeshFunctionGradFE< SCALAR_FE, SCALAR_COEFF > Class Template Reference

A MeshFunction representing the gradient of a function from a scalar finite element space (e.g. gradient of a solution of BVP). More...

#include <lf/fe/fe.h>

Public Types

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

Public Member Functions

 MeshFunctionGradFE (std::shared_ptr< const ScalarFESpace< SCALAR_FE > > fe_space, Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > dof_vector)
 Create a new MeshFunctionGradFE from a ScalarFESpace and a coefficient vector.
 
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 to which the original function belongs (i.e. before taking the gradient)
 
std::vector< Eigen::Matrix< Scalar, Eigen::Dynamic, 1 > > operator() (const lf::mesh::Entity &e, const Eigen::MatrixXd &local) const
 Evaluate the mesh function on a MeshEntity.
 

Private Attributes

std::shared_ptr< const ScalarFESpace< SCALAR_FE > > fe_space_
 Pointer to underlying finite element space.
 
Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > dof_vector_
 Reference to basis expansion coefficient vector for finite-element function
 

Detailed Description

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

A MeshFunction representing the gradient of a function from a scalar finite element space (e.g. gradient of a solution of BVP).

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

The MeshFunctionGradFE takes essentially two parameters:

Note
A MeshFunctionGradFE can be evaluated on entities of all codimensions except for points.
The gradient that is returned by this mesh function is w.r.t. to the global coordinate system of the mesh, i.e. the gradients of the shape function are multiplied by Geometry::JacobianInverseGramian(). Hence, what is returned is a "tangential gradient". For instance, if the evaluation operator is invoked for an edge of a 2D mesh, the gradient of the finite element function in the direction of that edge is returned. If the evaluation is requested for a (curved) triangle in 3D, then a gradient with 3 components is retrurned, which represents the tangential gradient of the finite element function in global 3D Euclidean coordinates.

Definition at line 47 of file mesh_function_grad_fe.h.

Member Typedef Documentation

◆ Scalar

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

Definition at line 49 of file mesh_function_grad_fe.h.

Constructor & Destructor Documentation

◆ MeshFunctionGradFE()

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

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

Parameters
fe_spacethe approximation space in which the function lies.
dof_vectorpasses the basis expansion coefficients
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_grad_fe.h.

Member Function Documentation

◆ getFESpace()

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

Convenience method to retrieve the finite element space to which the original function belongs (i.e. before taking the gradient)

Returns
shared_ptr to UniformScalarFESpace to which the original function belongs (i.e. before taking the gradient)

Definition at line 79 of file mesh_function_grad_fe.h.

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

◆ getMesh()

template<class SCALAR_FE , class SCALAR_COEFF >
std::shared_ptr< const mesh::Mesh > lf::fe::MeshFunctionGradFE< 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 69 of file mesh_function_grad_fe.h.

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

◆ operator()()

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

Evaluate the mesh function on a MeshEntity.

Parameters
ereference to a mesh entity
localmatrix whose columns contain reference coordinates of the evaluation points on the entity e
Returns
array of column vectors containing the gradients.

Definition at line 92 of file mesh_function_grad_fe.h.

References lf::base::RefEl::Dimension(), lf::fe::MeshFunctionGradFE< SCALAR_FE, SCALAR_COEFF >::dof_vector_, lf::fe::MeshFunctionGradFE< SCALAR_FE, SCALAR_COEFF >::fe_space_, lf::mesh::Entity::Geometry(), lf::geometry::Geometry::JacobianInverseGramian(), and lf::mesh::Entity::RefEl().

Member Data Documentation

◆ dof_vector_

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

Reference to basis expansion coefficient vector for finite-element function

Definition at line 125 of file mesh_function_grad_fe.h.

Referenced by lf::fe::MeshFunctionGradFE< SCALAR_FE, SCALAR_COEFF >::operator()().

◆ fe_space_

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

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