LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
segment_o2.h
1
9#ifndef SEGMENT_O2_H
10#define SEGMENT_O2_H
11
12#include "geometry_interface.h"
13
14namespace lf::geometry {
15
25class SegmentO2 : public Geometry {
26 public:
32 explicit SegmentO2(Eigen::Matrix<double, Eigen::Dynamic, 3> coords);
33
34 [[nodiscard]] dim_t DimLocal() const override { return 1; }
35 [[nodiscard]] dim_t DimGlobal() const override { return coords_.rows(); }
36 [[nodiscard]] base::RefEl RefEl() const override {
37 return base::RefEl::kSegment();
38 }
39
40 [[nodiscard]] Eigen::MatrixXd Global(
41 const Eigen::MatrixXd& local) const override;
42 [[nodiscard]] Eigen::MatrixXd Jacobian(
43 const Eigen::MatrixXd& local) const override;
44 [[nodiscard]] Eigen::MatrixXd JacobianInverseGramian(
45 const Eigen::MatrixXd& local) const override;
46 [[nodiscard]] Eigen::VectorXd IntegrationElement(
47 const Eigen::MatrixXd& local) const override;
48
49 [[nodiscard]] std::unique_ptr<Geometry> SubGeometry(dim_t codim,
50 dim_t i) const override;
51
52 [[nodiscard]] std::vector<std::unique_ptr<Geometry>> ChildGeometry(
53 const RefinementPattern& ref_pat, base::dim_t codim) const override;
54
55 private:
59 Eigen::Matrix<double, Eigen::Dynamic, 3> coords_;
60
61 /*
62 * SegmentO2 is parametrized by:
63 * alpha_ * x^2 + beta_ * x + gamma_
64 */
65 Eigen::Matrix<double, Eigen::Dynamic, 1> alpha_;
66 Eigen::Matrix<double, Eigen::Dynamic, 1> beta_;
67 Eigen::Matrix<double, Eigen::Dynamic, 1> gamma_;
68
69 /*
70 * Coefficients for efficient evaluation of JacobianInverseGramian() and
71 * IntegrationElement()
72 */
73 double alpha_squared_{0};
74 double alpha_beta_{0};
75 double beta_squared_{0};
76};
77
78} // namespace lf::geometry
79
80#endif /* SEGMENT_O2_H */
Represents a reference element with all its properties.
Definition ref_el.h:109
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition ref_el.h:153
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
A second-order segment in the plane or in 3D space.
Definition segment_o2.h:25
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Definition segment_o2.cc:66
SegmentO2(Eigen::Matrix< double, Eigen::Dynamic, 3 > coords)
Constructor building segment from vertex/midpoint coordinates.
Definition segment_o2.cc:15
std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pat, base::dim_t codim) const override
Generate geometry objects for child entities created in the course of refinement.
Definition segment_o2.cc:90
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition segment_o2.cc:28
Eigen::Matrix< double, Eigen::Dynamic, 3 > coords_
Coordinates of the 3 vertices/midpoints, stored in matrix columns.
Definition segment_o2.h:59
Eigen::Matrix< double, Eigen::Dynamic, 1 > beta_
Definition segment_o2.h:66
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition segment_o2.cc:46
dim_t DimLocal() const override
Dimension of the domain of this mapping.
Definition segment_o2.h:34
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Definition segment_o2.h:36
Eigen::Matrix< double, Eigen::Dynamic, 1 > gamma_
Definition segment_o2.h:67
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition segment_o2.h:35
Eigen::Matrix< double, Eigen::Dynamic, 1 > alpha_
Definition segment_o2.h:65
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 segment_o2.cc:78
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition segment_o2.cc:39
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