LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
tria_o1.h
1#ifndef INCGb84bb391fa744cdd9159293bfcb5b311
2#define INCGb84bb391fa744cdd9159293bfcb5b311
3
4#include "geometry_interface.h"
5
6namespace lf::geometry {
7
17 const Eigen::Matrix<double, Eigen::Dynamic, 3>& coords,
18 double tol = 1.0E-8);
19
23class TriaO1 : public Geometry {
24 public:
25 explicit TriaO1(Eigen::Matrix<double, Eigen::Dynamic, 3> coords);
26
27 [[nodiscard]] dim_t DimLocal() const override { return 2; }
28 [[nodiscard]] dim_t DimGlobal() const override { return coords_.rows(); }
29 [[nodiscard]] base::RefEl RefEl() const override {
30 return base::RefEl::kTria();
31 }
32 [[nodiscard]] Eigen::MatrixXd Global(
33 const Eigen::MatrixXd& local) const override;
34
35 [[nodiscard]] Eigen::MatrixXd Jacobian(
36 const Eigen::MatrixXd& local) const override {
37 return jacobian_.replicate(1, local.cols());
38 }
39 [[nodiscard]] Eigen::MatrixXd JacobianInverseGramian(
40 const Eigen::MatrixXd& local) const override {
41 return jacobian_inverse_gramian_.replicate(1, local.cols());
42 }
43 [[nodiscard]] Eigen::VectorXd IntegrationElement(
44 const Eigen::MatrixXd& local) const override {
45 return Eigen::VectorXd::Constant(local.cols(), integrationElement_);
46 }
47 [[nodiscard]] std::unique_ptr<Geometry> SubGeometry(dim_t codim,
48 dim_t i) const override;
49
58 [[nodiscard]] std::vector<std::unique_ptr<Geometry>> ChildGeometry(
59 const RefinementPattern& ref_pat, lf::base::dim_t codim) const override;
60
61 private:
62 Eigen::Matrix<double, Eigen::Dynamic, 3> coords_;
63 Eigen::Matrix<double, Eigen::Dynamic, 2> jacobian_;
64 Eigen::Matrix<double, Eigen::Dynamic, 2> jacobian_inverse_gramian_;
66};
67} // namespace lf::geometry
68
69#endif // INCGb84bb391fa744cdd9159293bfcb5b311
Represents a reference element with all its properties.
Definition ref_el.h:109
static constexpr RefEl kTria()
Returns the reference triangle.
Definition ref_el.h:161
Interface class for shape information on a mesh cell in the spirit of parametric finite element metho...
base::RefEl::dim_t dim_t
Abstract interface class for encoding topological local refinement
An affine triangle in the plane or in 3D space.
Definition tria_o1.h:23
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition tria_o1.h:35
TriaO1(Eigen::Matrix< double, Eigen::Dynamic, 3 > coords)
Definition tria_o1.cc:50
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition tria_o1.cc:70
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition tria_o1.h:28
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition tria_o1.h:39
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_
Definition tria_o1.h:63
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Definition tria_o1.h:43
std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pat, lf::base::dim_t codim) const override
creation of child geometries as specified in refinement pattern
Definition tria_o1.cc:99
dim_t DimLocal() const override
Dimension of the domain of this mapping.
Definition tria_o1.h:27
Eigen::Matrix< double, Eigen::Dynamic, 3 > coords_
Definition tria_o1.h:62
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 tria_o1.cc:76
double integrationElement_
Definition tria_o1.h:65
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_inverse_gramian_
Definition tria_o1.h:64
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Definition tria_o1.h:29
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 assertNonDegenerateTriangle(const Eigen::Matrix< double, Eigen::Dynamic, 3 > &coords, double tol)
Asserting non-degenerate shape of a flat triangle.
Definition tria_o1.cc:10