LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
quad_o1.h
1#ifndef INCGb193e889f4e245959d124b0ded8a3b68
2#define INCGb193e889f4e245959d124b0ded8a3b68
3
4#include "geometry_interface.h"
5
6namespace lf::geometry {
7
17 const Eigen::Matrix<double, Eigen::Dynamic, 4>& coords,
18 double tol = 1.0E-8);
19
28class QuadO1 : public Geometry {
29 public:
38 explicit QuadO1(Eigen::Matrix<double, Eigen::Dynamic, 4> coords);
39
40 [[nodiscard]] dim_t DimLocal() const override { return 2; }
41 [[nodiscard]] dim_t DimGlobal() const override { return coords_.rows(); }
42 [[nodiscard]] base::RefEl RefEl() const override {
43 return base::RefEl::kQuad();
44 }
45
50 [[nodiscard]] Eigen::MatrixXd Global(
51 const Eigen::MatrixXd& local) const override;
56 [[nodiscard]] Eigen::MatrixXd Jacobian(
57 const Eigen::MatrixXd& local) const override;
62 [[nodiscard]] Eigen::MatrixXd JacobianInverseGramian(
63 const Eigen::MatrixXd& local) const override;
68 [[nodiscard]] Eigen::VectorXd IntegrationElement(
69 const Eigen::MatrixXd& local) const override;
70
72 [[nodiscard]] std::unique_ptr<Geometry> SubGeometry(dim_t codim,
73 dim_t i) const override;
74
76 [[nodiscard]] bool isAffine() const override { return false; }
77
84 [[nodiscard]] std::vector<std::unique_ptr<Geometry>> ChildGeometry(
85 const RefinementPattern& ref_pat, lf::base::dim_t codim) const override;
86
87 private:
89 Eigen::Matrix<double, Eigen::Dynamic, 4> coords_;
90};
91
97class Parallelogram : public Geometry {
98 public:
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);
121
122 [[nodiscard]] dim_t DimLocal() const override { return 2; }
123 [[nodiscard]] dim_t DimGlobal() const override { return coords_.rows(); }
124 [[nodiscard]] base::RefEl RefEl() const override {
125 return base::RefEl::kQuad();
126 }
127
128 [[nodiscard]] Eigen::MatrixXd Global(
129 const Eigen::MatrixXd& local) const override;
130 [[nodiscard]] Eigen::MatrixXd Jacobian(
131 const Eigen::MatrixXd& local) const override;
132 [[nodiscard]] Eigen::MatrixXd JacobianInverseGramian(
133 const Eigen::MatrixXd& local) const override;
134 [[nodiscard]] Eigen::VectorXd IntegrationElement(
135 const Eigen::MatrixXd& local) const override;
136
138 [[nodiscard]] std::unique_ptr<Geometry> SubGeometry(dim_t codim,
139 dim_t i) const override;
140
146 [[nodiscard]] bool isAffine() const override { return true; }
147
154 [[nodiscard]] std::vector<std::unique_ptr<Geometry>> ChildGeometry(
155 const RefinementPattern& ref_pat, lf::base::dim_t codim) const override;
156
157 private:
159 void init();
160
167 Eigen::Matrix<double, Eigen::Dynamic, 4> coords_;
170 Eigen::Matrix<double, Eigen::Dynamic, 2> jacobian_;
172 Eigen::Matrix<double, Eigen::Dynamic, 2> jacobian_inverse_gramian_;
175};
176
177} // namespace lf::geometry
178
179#endif // INCGb193e889f4e245959d124b0ded8a3b68
Represents a reference element with all its properties.
Definition ref_el.h:109
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition ref_el.h:169
Interface class for shape information on a mesh cell in the spirit of parametric finite element metho...
base::RefEl::dim_t dim_t
Affine quadrilateral = parallelogram.
Definition quad_o1.h:97
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 quad_o1.cc:327
double integrationElement_
constant Gramian determinant
Definition quad_o1.h:174
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_inverse_gramian_
Constant matrix ( )
Definition quad_o1.h:172
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition quad_o1.cc:306
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition quad_o1.cc:316
dim_t DimLocal() const override
Dimension of the domain of this mapping.
Definition quad_o1.h:122
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_
Matrix of affine mapping generating the parallelogram, agrees with constant Jacobian.
Definition quad_o1.h:170
bool isAffine() const override
element shape by affine mapping from reference element
Definition quad_o1.h:146
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition quad_o1.h:123
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Definition quad_o1.cc:321
Eigen::Matrix< double, Eigen::Dynamic, 4 > coords_
Coordinates of the a four vertices, stored in matrix columns.
Definition quad_o1.h:167
Parallelogram(Eigen::Matrix< double, Eigen::Dynamic, 4 > coords)
Constructor building a parallelogram from four vertex coordinates.
Definition quad_o1.cc:258
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition quad_o1.cc:312
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Definition quad_o1.h:124
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.
Definition quad_o1.cc:351
void init()
performs initialization of data members
Definition quad_o1.cc:290
Bilinear quadrilateral element shape.
Definition quad_o1.h:28
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition quad_o1.cc:67
bool isAffine() const override
element shape by affine mapping from reference element
Definition quad_o1.h:76
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 quad_o1.cc:159
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Definition quad_o1.cc:136
QuadO1(Eigen::Matrix< double, Eigen::Dynamic, 4 > coords)
Constructor building quadrilateral from vertex coordinates.
Definition quad_o1.cc:60
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition quad_o1.cc:106
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition quad_o1.cc:85
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition quad_o1.h:41
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.
Definition quad_o1.cc:181
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Definition quad_o1.h:42
Eigen::Matrix< double, Eigen::Dynamic, 4 > coords_
Coordinates of the a four vertices, stored in matrix columns.
Definition quad_o1.h:89
dim_t DimLocal() const override
Dimension of the domain of this mapping.
Definition quad_o1.h:40
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
bool assertNonDegenerateQuad(const Eigen::Matrix< double, Eigen::Dynamic, 4 > &coords, double tol)
Asserting a non-degenerate bilinear quadrilateral.
Definition quad_o1.cc:12