LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
point.h
1#ifndef INCG184598b89ca44fe1a1e7a043bc32da06
2#define INCG184598b89ca44fe1a1e7a043bc32da06
3
4#include <lf/base/base.h>
5
6#include "geometry_interface.h"
7
8namespace lf::geometry {
9
10class Point : public Geometry {
11 public:
12 explicit Point(Eigen::VectorXd coord) : coord_(std::move(coord)) {}
13
14 [[nodiscard]] dim_t DimLocal() const override { return 0; }
15
16 [[nodiscard]] dim_t DimGlobal() const override { return coord_.rows(); }
17
18 [[nodiscard]] base::RefEl RefEl() const override {
19 return base::RefEl::kPoint();
20 }
21
22 [[nodiscard]] Eigen::MatrixXd Global(
23 const Eigen::MatrixXd& local) const override;
24
25 [[nodiscard]] Eigen::MatrixXd Jacobian(
26 const Eigen::MatrixXd& local) const override;
27 [[nodiscard]] Eigen::MatrixXd JacobianInverseGramian(
28 const Eigen::MatrixXd& local) const override;
29 [[nodiscard]] Eigen::VectorXd IntegrationElement(
30 const Eigen::MatrixXd& local) const override;
31 [[nodiscard]] std::unique_ptr<Geometry> SubGeometry(dim_t codim,
32 dim_t i) const override;
33
37 [[nodiscard]] std::vector<std::unique_ptr<Geometry>> ChildGeometry(
38 const RefinementPattern& ref_pattern, base::dim_t codim) const override;
39
40 private:
41 Eigen::VectorXd coord_;
42};
43
44} // namespace lf::geometry
45
46#endif // INCG184598b89ca44fe1a1e7a043bc32da06
Represents a reference element with all its properties.
Definition ref_el.h:109
static constexpr RefEl kPoint()
Returns the (0-dimensional) reference point.
Definition ref_el.h:144
Interface class for shape information on a mesh cell in the spirit of parametric finite element metho...
base::RefEl::dim_t dim_t
std::unique_ptr< Geometry > SubGeometry(dim_t codim, dim_t i) const override
Construct a new Geometry() object that describes the geometry of the i-th sub-entity with codimension...
Definition point.cc:23
Point(Eigen::VectorXd coord)
Definition point.h:12
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition point.cc:10
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition point.cc:4
Eigen::VectorXd coord_
Definition point.h:41
std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pattern, base::dim_t codim) const override
the child geometry is just a copy of the point geometry
Definition point.cc:30
dim_t DimLocal() const override
Dimension of the domain of this mapping.
Definition point.h:14
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Definition point.h:18
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition point.h:16
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition point.cc:14
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Definition point.cc:19
Abstract interface class for encoding topological local refinement
unsigned int dim_t
type for dimensions and co-dimensions and numbers derived from them
Definition types.h:32
Defines the Geometry interface and provides a number of classes that implement this interface + addit...
Definition compose.h:13