LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
mass_matrix_provider.cc
1#include "mass_matrix_provider.h"
2
3#include <lf/uscalfe/lagr_fe.h>
4
5#include <cmath>
6#include <iomanip>
7#include <iostream>
8
10
11Eigen::MatrixXd MassMatrixProvider::Eval(const lf::mesh::Entity &entity) const {
12 // Only triangles are supported
13 LF_VERIFY_MSG(entity.RefEl() == lf::base::RefEl::kTria(),
14 "Unsupported cell type " << entity.RefEl());
15
16 // Get the geometry of the entity
17 const auto *geom = entity.Geometry();
18
19 // Compute the element matrix for the product of baricentric functions
20 Eigen::MatrixXd elem_mat(3, 3);
21 // clang-format off
22 elem_mat << 2, 1, 1,
23 1, 2, 1,
24 1, 1, 2;
25 // clang-format on
26 elem_mat *= lf::geometry::Volume(*geom) / 12.0;
27
28 return elem_mat;
29}
30
31} // namespace projects::hldo_sphere::assemble
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 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.
Eigen::MatrixXd Eval(const lf::mesh::Entity &entity) const
Compute the element matrix for a given cell of a mesh.
double Volume(const Geometry &geo)
Compute the (approximate) volume (area) of a shape.
Collection of matrix and vector element providers.
Definition: assemble.h:26