LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
segment_o1.h
1#ifndef INCG0ec30f6fc5554ec9bea4d9b83e83ea42
2#define INCG0ec30f6fc5554ec9bea4d9b83e83ea42
3
4#include "geometry_interface.h"
5
6namespace lf::geometry {
7
11class SegmentO1 : public Geometry {
12 public:
13 explicit SegmentO1(Eigen::Matrix<double, Eigen::Dynamic, 2> coords)
14 : coords_(std::move(coords)) {}
15
16 [[nodiscard]] dim_t DimLocal() const override { return 1; }
17 [[nodiscard]] dim_t DimGlobal() const override { return coords_.rows(); }
18 [[nodiscard]] base::RefEl RefEl() const override {
19 return base::RefEl::kSegment();
20 }
21 [[nodiscard]] Eigen::MatrixXd Global(
22 const Eigen::MatrixXd& local) const override;
23 [[nodiscard]] Eigen::MatrixXd Jacobian(
24 const Eigen::MatrixXd& local) const override;
25 [[nodiscard]] Eigen::MatrixXd JacobianInverseGramian(
26 const Eigen::MatrixXd& local) const override;
27 [[nodiscard]] Eigen::VectorXd IntegrationElement(
28 const Eigen::MatrixXd& local) const override;
29 [[nodiscard]] std::unique_ptr<Geometry> SubGeometry(dim_t codim,
30 dim_t i) const override;
31 [[nodiscard]] std::vector<std::unique_ptr<Geometry>> ChildGeometry(
32 const RefinementPattern& ref_pat, base::dim_t codim) const override;
33
34 private:
35 Eigen::Matrix<double, Eigen::Dynamic, 2> coords_;
36};
37
38} // namespace lf::geometry
39
40#endif // INCG0ec30f6fc5554ec9bea4d9b83e83ea42
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 straight edge defined by the location of its two endpoints.
Definition segment_o1.h:11
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition segment_o1.cc:11
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition segment_o1.cc:15
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_o1.cc:33
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Definition segment_o1.h:18
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_o1.cc:27
SegmentO1(Eigen::Matrix< double, Eigen::Dynamic, 2 > coords)
Definition segment_o1.h:13
Eigen::Matrix< double, Eigen::Dynamic, 2 > coords_
Definition segment_o1.h:35
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition segment_o1.h:17
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition segment_o1.cc:7
dim_t DimLocal() const override
Dimension of the domain of this mapping.
Definition segment_o1.h:16
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_o1.cc:45
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