12#include "geometry_interface.h"
32 explicit SegmentO2(Eigen::Matrix<double, Eigen::Dynamic, 3> coords);
40 [[nodiscard]] Eigen::MatrixXd
Global(
41 const Eigen::MatrixXd& local)
const override;
42 [[nodiscard]] Eigen::MatrixXd
Jacobian(
43 const Eigen::MatrixXd& local)
const override;
45 const Eigen::MatrixXd& local)
const override;
47 const Eigen::MatrixXd& local)
const override;
49 [[nodiscard]] std::unique_ptr<Geometry>
SubGeometry(dim_t codim,
50 dim_t i)
const override;
52 [[nodiscard]] std::vector<std::unique_ptr<Geometry>>
ChildGeometry(
59 Eigen::Matrix<double, Eigen::Dynamic, 3>
coords_;
65 Eigen::Matrix<double, Eigen::Dynamic, 1>
alpha_;
66 Eigen::Matrix<double, Eigen::Dynamic, 1>
beta_;
67 Eigen::Matrix<double, Eigen::Dynamic, 1>
gamma_;
Represents a reference element with all its properties.
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Interface class for shape information on a mesh cell in the spirit of parametric finite element metho...
Abstract interface class for encoding topological local refinement
A second-order segment in the plane or in 3D space.
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
SegmentO2(Eigen::Matrix< double, Eigen::Dynamic, 3 > coords)
Constructor building segment from vertex/midpoint coordinates.
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.
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Eigen::Matrix< double, Eigen::Dynamic, 3 > coords_
Coordinates of the 3 vertices/midpoints, stored in matrix columns.
Eigen::Matrix< double, Eigen::Dynamic, 1 > beta_
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.
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Eigen::Matrix< double, Eigen::Dynamic, 1 > gamma_
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Eigen::Matrix< double, Eigen::Dynamic, 1 > alpha_
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::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
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...