LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
mesh_function_whitney_zero.h
1#ifndef HLDO_SPHERE_MESH_FUNCTION_WHITNEY_ZERO_H
2#define HLDO_SPHERE_MESH_FUNCTION_WHITNEY_ZERO_H
3/*
4 *@file mesh_function_whitney_zero.h
5 */
6#include <lf/uscalfe/uscalfe.h>
7
9
18template <typename SCALAR>
20 public:
35 MeshFunctionWhitneyZero(const Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>& mu,
36 const std::shared_ptr<const lf::mesh::Mesh> mesh)
37 : mu_(mu), mesh_(mesh) {
38 // make sure mu has the right size
39 const lf::base::size_type n = mesh->NumEntities(2);
40 const int mu_size = mu.size();
41 LF_VERIFY_MSG(n == mu_size,
42 "Not the right number of basis expansion coefficiants in "
43 "mu expected: "
44 << n << " given: " << mu_size);
45 };
46
75 std::vector<SCALAR> operator()(const lf::mesh::Entity& e,
76 const Eigen::MatrixXd& local) const {
77 // Only triangles are supported
78 LF_VERIFY_MSG(e.RefEl() == lf::base::RefEl::kTria(),
79 "Unsupported cell type " << e.RefEl());
80
81 // get the lambdas
83
84 // get lambda values
85 const Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> lambda =
86 hat_func.EvalReferenceShapeFunctions(local);
87
88 // define return vector
89 std::vector<SCALAR> vals(local.cols());
90
91 // get global indices of the edges
92 std::vector<lf::base::size_type> global_indices(3);
93 auto points = e.SubEntities(2);
94 for (lf::base::size_type i = 0; i < 3; i++) {
95 global_indices[i] = mesh_->Index(*points[i]);
96 }
97
98 // for each point evaluate the basis functions and add the
99 // values together
100 for (int i = 0; i < local.cols(); i++) {
101 vals[i] = mu_(global_indices[0]) * lambda(0) +
102 mu_(global_indices[1]) * lambda(1) +
103 mu_(global_indices[2]) * lambda(2);
104 }
105
106 return vals;
107 }
108
109 private:
110 const Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>& mu_;
111 const std::shared_ptr<const lf::mesh::Mesh> mesh_;
112};
113
114} // namespace projects::hldo_sphere::post_processing
115#endif // HLDO_SPHERE_MESH_FUNCTION_WHITNEY_ONE_H
static constexpr RefEl kTria()
Returns the reference triangle.
Definition: ref_el.h:158
Interface class representing a topological entity in a cellular complex
Definition: entity.h:39
virtual nonstd::span< const Entity *const > SubEntities(unsigned rel_codim) const =0
Return all sub entities of this entity that have the given codimension (w.r.t. this entity!...
virtual base::RefEl RefEl() const =0
Describes the reference element type of this entity.
Linear Lagrange finite element on triangular reference element.
Definition: lagr_fe.h:58
Eigen::Matrix< SCALAR, Eigen::Dynamic, Eigen::Dynamic > EvalReferenceShapeFunctions(const Eigen::MatrixXd &refcoords) const override
Evaluation of all reference shape functions in a number of points.
Definition: lagr_fe.h:104
Provides Mesh Function for basis expansion coefficients of the zero form.
MeshFunctionWhitneyZero(const Eigen::Matrix< SCALAR, Eigen::Dynamic, 1 > &mu, const std::shared_ptr< const lf::mesh::Mesh > mesh)
basic constructor
unsigned int size_type
general type for variables related to size of arrays
Definition: base.h:24
Postprocessing such as computing norms and Mesh Functions and plot scripts.