LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
geometry_interface.h
1#ifndef INCG7ed6b0d4d9244155819c464fc4eb9bbb
2#define INCG7ed6b0d4d9244155819c464fc4eb9bbb
3
4#include <lf/base/base.h>
5
6#include <memory>
7
8#include "refinement_pattern.h"
9
10namespace lf::geometry {
11
21class Geometry {
22 protected:
23 Geometry() = default;
24 Geometry(const Geometry&) = default;
25 Geometry(Geometry&&) = default;
26 Geometry& operator=(const Geometry&) = default;
28
29 public:
31
35 [[nodiscard]] virtual dim_t DimLocal() const = 0;
36
40 [[nodiscard]] virtual dim_t DimGlobal() const = 0;
41
45 [[nodiscard]] virtual base::RefEl RefEl() const = 0;
46
82 [[nodiscard]] virtual Eigen::MatrixXd Global(
83 const Eigen::MatrixXd& local) const = 0;
84
101 [[nodiscard]] virtual Eigen::MatrixXd Jacobian(
102 const Eigen::MatrixXd& local) const = 0;
103
131 [[nodiscard]] virtual Eigen::MatrixXd JacobianInverseGramian(
132 const Eigen::MatrixXd& local) const = 0;
133
157 [[nodiscard]] virtual Eigen::VectorXd IntegrationElement(
158 const Eigen::MatrixXd& local) const = 0;
159
180 [[nodiscard]] virtual std::unique_ptr<Geometry> SubGeometry(
181 dim_t codim, dim_t i) const = 0;
182
204 [[nodiscard]] virtual std::vector<std::unique_ptr<Geometry>> ChildGeometry(
205 const RefinementPattern& ref_pat, lf::base::dim_t codim) const = 0;
206
220 [[nodiscard]] virtual bool isAffine() const { return true; }
221
225 virtual ~Geometry() = default;
226
227}; // class Geometry
228
237double Volume(const Geometry& geo);
238
261inline Eigen::MatrixXd Corners(const Geometry& geo) {
262 return geo.Global(geo.RefEl().NodeCoords());
263}
264
265} // namespace lf::geometry
266
267#endif // INCG7ed6b0d4d9244155819c464fc4eb9bbb
Represents a reference element with all its properties.
Definition ref_el.h:109
unsigned int dim_t
Definition ref_el.h:132
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
Definition ref_el.h:241
Interface class for shape information on a mesh cell in the spirit of parametric finite element metho...
virtual base::RefEl RefEl() const =0
The Reference element that defines the domain of this mapping.
virtual std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pat, lf::base::dim_t codim) const =0
Generate geometry objects for child entities created in the course of refinement.
base::RefEl::dim_t dim_t
virtual Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const =0
Map a number of points in local coordinates into the global coordinate system.
virtual dim_t DimLocal() const =0
Dimension of the domain of this mapping.
virtual dim_t DimGlobal() const =0
Dimension of the image of this mapping.
Geometry & operator=(const Geometry &)=default
virtual Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const =0
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
virtual Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const =0
The integration element (factor appearing in integral transformation formula, see below) at number of...
virtual std::unique_ptr< Geometry > SubGeometry(dim_t codim, dim_t i) const =0
Construct a new Geometry() object that describes the geometry of the i-th sub-entity with codimension...
virtual bool isAffine() const
element shape by affine mapping from reference element
virtual ~Geometry()=default
Virtual destructor.
Geometry(const Geometry &)=default
virtual Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const =0
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Geometry & operator=(Geometry &&)=default
Geometry(Geometry &&)=default
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
Eigen::MatrixXd Corners(const Geometry &geo)
The corners of a shape with piecewise smooth boundary.
double Volume(const Geometry &geo)
Compute the (approximate) volume (area) of a shape.