LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
mesh_function_grad_fe.h
1
9#ifndef __b6997524e2834b5b8e4bba019fb35cc6
10#define __b6997524e2834b5b8e4bba019fb35cc6
11#include "scalar_fe_space.h"
12
13namespace lf::fe {
14
45template <class SCALAR_FE, class SCALAR_COEFF>
47 public:
48 using Scalar = decltype(SCALAR_FE(0) * SCALAR_COEFF(0));
49
60 MeshFunctionGradFE(std::shared_ptr<const ScalarFESpace<SCALAR_FE>> fe_space,
61 Eigen::Matrix<SCALAR_COEFF, Eigen::Dynamic, 1> dof_vector)
62 : fe_space_(std::move(fe_space)), dof_vector_(std::move(dof_vector)) {}
63
68 [[nodiscard]] std::shared_ptr<const mesh::Mesh> getMesh() const {
69 return fe_space_->Mesh();
70 }
71
78 [[nodiscard]] std::shared_ptr<const ScalarFESpace<SCALAR_FE>> getFESpace()
79 const {
80 return fe_space_;
81 }
82
91 std::vector<Eigen::Matrix<Scalar, Eigen::Dynamic, 1>> operator()(
92 const lf::mesh::Entity& e, const Eigen::MatrixXd& local) const {
93 // Access information on local shape functions for the entity
94 auto grad_sf_eval =
95 fe_space_->ShapeFunctionLayout(e)->GradientsReferenceShapeFunctions(
96 local);
97 // Fetch the coefficients of the global shape functions associated with
98 // the current mesh entity
99 Eigen::Matrix<SCALAR_COEFF, 1, Eigen::Dynamic> local_dofs(
100 1, grad_sf_eval.rows());
101 auto global_dofs = fe_space_->LocGlobMap().GlobalDofIndices(e);
102 for (Eigen::Index i = 0; i < grad_sf_eval.rows(); ++i) {
103 local_dofs(i) = dof_vector_(global_dofs[i]);
104 }
105 // gradients w.r.t. reference element coordinates \hat{x}
106 auto local_grads = (local_dofs * grad_sf_eval).eval();
107 // Transform to Cartesian coordinates
108 auto jac_t = e.Geometry()->JacobianInverseGramian(local);
109 auto dim_local = e.RefEl().Dimension();
110 std::vector<Eigen::Matrix<Scalar, Eigen::Dynamic, 1>> result(local.cols());
111 // Transform all the local gradients in the evaluation points
112 for (std::size_t i = 0; i < result.size(); ++i) {
113 result[i] = jac_t.block(0, dim_local * i, jac_t.rows(), dim_local) *
114 local_grads.block(0, i * dim_local, 1, dim_local).transpose();
115 }
116 return result;
117 }
118
119 private:
121 std::shared_ptr<const ScalarFESpace<SCALAR_FE>> fe_space_;
124 Eigen::Matrix<SCALAR_COEFF, Eigen::Dynamic, 1> dof_vector_;
125};
126
127// deduction guide
128template <class T, class SCALAR_COEFF>
129MeshFunctionGradFE(std::shared_ptr<T>,
130 const Eigen::Matrix<SCALAR_COEFF, Eigen::Dynamic, 1>&)
132
133} // namespace lf::fe
134
135#endif // __b6997524e2834b5b8e4bba019fb35cc6
constexpr dim_t Dimension() const
Return the dimension of this reference element.
Definition: ref_el.h:201
A MeshFunction representing the gradient of a function from a finite element space (e....
std::shared_ptr< const ScalarFESpace< SCALAR_FE > > fe_space_
Pointer to underlying finite element space.
std::shared_ptr< const mesh::Mesh > getMesh() 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....
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.
Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > dof_vector_
Reference to basis expansion coefficient vector for finite-element function
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 ScalarUniformFESpace and a coefficient vector.
decltype(SCALAR_FE(0) *SCALAR_COEFF(0)) Scalar
Space of scalar valued finite element functions on a hybrid 2D mesh
virtual Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const =0
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Interface class representing a topological entity in a cellular complex
Definition: entity.h:39
virtual const geometry::Geometry * Geometry() const =0
Describes the geometry of this entity.
virtual base::RefEl RefEl() const =0
Describes the reference element type of this entity.
Collects data structures and algorithms designed for scalar finite element methods primarily meant fo...
Definition: fe.h:47
MeshFunctionGradFE(std::shared_ptr< T >, const Eigen::Matrix< SCALAR_COEFF, Eigen::Dynamic, 1 > &) -> MeshFunctionGradFE< typename T::Scalar, SCALAR_COEFF >