1#ifndef INCG7ed6b0d4d9244155819c464fc4eb9bbb
2#define INCG7ed6b0d4d9244155819c464fc4eb9bbb
4#include <lf/base/base.h>
8#include "refinement_pattern.h"
82 [[nodiscard]]
virtual Eigen::MatrixXd
Global(
83 const Eigen::MatrixXd& local)
const = 0;
102 const Eigen::MatrixXd& local)
const = 0;
132 const Eigen::MatrixXd& local)
const = 0;
158 const Eigen::MatrixXd& local)
const = 0;
220 [[nodiscard]]
virtual bool isAffine()
const {
return true; }
Represents a reference element with all its properties.
const Eigen::MatrixXd & NodeCoords() const
Get the coordinates of the nodes of this reference element.
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.
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
Defines the Geometry interface and provides a number of classes that implement this interface + addit...
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.