1#ifndef INCGb193e889f4e245959d124b0ded8a3b68
2#define INCGb193e889f4e245959d124b0ded8a3b68
4#include "geometry_interface.h"
17 const Eigen::Matrix<double, Eigen::Dynamic, 4>& coords,
38 explicit QuadO1(Eigen::Matrix<double, Eigen::Dynamic, 4> coords);
50 [[nodiscard]] Eigen::MatrixXd
Global(
51 const Eigen::MatrixXd& local)
const override;
56 [[nodiscard]] Eigen::MatrixXd
Jacobian(
57 const Eigen::MatrixXd& local)
const override;
63 const Eigen::MatrixXd& local)
const override;
69 const Eigen::MatrixXd& local)
const override;
72 [[nodiscard]] std::unique_ptr<Geometry>
SubGeometry(dim_t codim,
73 dim_t i)
const override;
76 [[nodiscard]]
bool isAffine()
const override {
return false; }
84 [[nodiscard]] std::vector<std::unique_ptr<Geometry>>
ChildGeometry(
89 Eigen::Matrix<double, Eigen::Dynamic, 4>
coords_;
108 explicit Parallelogram(Eigen::Matrix<double, Eigen::Dynamic, 4> coords);
119 explicit Parallelogram(
const Eigen::VectorXd& p0,
const Eigen::VectorXd& p1,
120 const Eigen::VectorXd& p2);
128 [[nodiscard]] Eigen::MatrixXd
Global(
129 const Eigen::MatrixXd& local)
const override;
130 [[nodiscard]] Eigen::MatrixXd
Jacobian(
131 const Eigen::MatrixXd& local)
const override;
133 const Eigen::MatrixXd& local)
const override;
135 const Eigen::MatrixXd& local)
const override;
138 [[nodiscard]] std::unique_ptr<Geometry>
SubGeometry(dim_t codim,
139 dim_t i)
const override;
146 [[nodiscard]]
bool isAffine()
const override {
return true; }
154 [[nodiscard]] std::vector<std::unique_ptr<Geometry>>
ChildGeometry(
167 Eigen::Matrix<double, Eigen::Dynamic, 4>
coords_;
Represents a reference element with all its properties.
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Interface class for shape information on a mesh cell in the spirit of parametric finite element metho...
Affine quadrilateral = parallelogram.
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...
double integrationElement_
constant Gramian determinant
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_inverse_gramian_
Constant matrix ( )
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
dim_t DimLocal() const override
Dimension of the domain of this mapping.
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_
Matrix of affine mapping generating the parallelogram, agrees with constant Jacobian.
bool isAffine() const override
element shape by affine mapping from reference element
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Eigen::Matrix< double, Eigen::Dynamic, 4 > coords_
Coordinates of the a four vertices, stored in matrix columns.
Parallelogram(Eigen::Matrix< double, Eigen::Dynamic, 4 > coords)
Constructor building a parallelogram from four vertex coordinates.
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pat, lf::base::dim_t codim) const override
Generate geometry objects for child entities created in the course of refinement.
void init()
performs initialization of data members
Bilinear quadrilateral element shape.
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
bool isAffine() const override
element shape by affine mapping from reference element
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...
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
QuadO1(Eigen::Matrix< double, Eigen::Dynamic, 4 > coords)
Constructor building quadrilateral from vertex coordinates.
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
dim_t DimGlobal() const override
Dimension of the image of this mapping.
std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pat, lf::base::dim_t codim) const override
Generate geometry objects for child entities created in the course of refinement.
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Eigen::Matrix< double, Eigen::Dynamic, 4 > coords_
Coordinates of the a four vertices, stored in matrix columns.
dim_t DimLocal() const override
Dimension of the domain of this mapping.
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...
bool assertNonDegenerateQuad(const Eigen::Matrix< double, Eigen::Dynamic, 4 > &coords, double tol)
Asserting a non-degenerate bilinear quadrilateral.