1#ifndef HLDO_SPHERE_MESH_FUNCTION_WHITNEY_ZERO_H
2#define HLDO_SPHERE_MESH_FUNCTION_WHITNEY_ZERO_H
6#include <lf/uscalfe/uscalfe.h>
18template <
typename SCALAR>
36 const std::shared_ptr<const lf::mesh::Mesh> mesh)
37 : mu_(mu), mesh_(mesh) {
40 const int mu_size = mu.size();
41 LF_VERIFY_MSG(n == mu_size,
42 "Not the right number of basis expansion coefficiants in "
44 << n <<
" given: " << mu_size);
76 const Eigen::MatrixXd& local)
const {
79 "Unsupported cell type " << e.
RefEl());
85 const Eigen::Matrix<SCALAR, Eigen::Dynamic, Eigen::Dynamic> lambda =
89 std::vector<SCALAR> vals(local.cols());
92 std::vector<lf::base::size_type> global_indices(3);
95 global_indices[i] = mesh_->Index(*points[i]);
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);
110 const Eigen::Matrix<SCALAR, Eigen::Dynamic, 1>& mu_;
111 const std::shared_ptr<const lf::mesh::Mesh> mesh_;
static constexpr RefEl kTria()
Returns the reference triangle.
Interface class representing a topological entity in a cellular complex
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.
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.
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
Postprocessing such as computing norms and Mesh Functions and plot scripts.