LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
segment_o2.cc
1
9#include "segment_o2.h"
10
11#include "point.h"
12
13namespace lf::geometry {
14
15SegmentO2::SegmentO2(Eigen::Matrix<double, Eigen::Dynamic, 3> coords)
16 : coords_(std::move(coords)) {
17 // polynomial of degree 2: alpha * x^2 + beta * x + gamma
18 alpha_ = 2. * (coords_.col(1) + coords_.col(0)) - 4. * coords_.col(2);
19 beta_ = 4. * coords_.col(2) - 3. * coords_.col(0) - coords_.col(1);
20 gamma_ = coords_.col(0);
21
22 // coefficients for JacobianInverseGramian and IntegrationElement
23 alpha_squared_ = alpha_.squaredNorm();
25 beta_squared_ = beta_.squaredNorm();
26}
27
28Eigen::MatrixXd SegmentO2::Global(const Eigen::MatrixXd& local) const {
29 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
30 "local coordinates out of bounds for reference element");
31
32 // evaluate polynomial with Horner scheme
33 Eigen::VectorXd local_vec = local.transpose();
34 auto tmp = ((alpha_ * local).colwise() + beta_).array().rowwise();
35
36 return (tmp * local_vec.transpose().array()).matrix().colwise() + gamma_;
37}
38
39Eigen::MatrixXd SegmentO2::Jacobian(const Eigen::MatrixXd& local) const {
40 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
41 "local coordinates out of bounds for reference element");
42
43 return (2. * alpha_ * local).colwise() + beta_;
44}
45
47 const Eigen::MatrixXd& local) const {
48 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
49 "local coordinates out of bounds for reference element");
50
51 auto jacobian = (2. * alpha_ * local).colwise() + beta_;
52
53 if (DimGlobal() == 1) {
54 return jacobian.cwiseInverse();
55 }
56
57 // evaluate polynomial with Horner scheme
58 const auto jTj =
59 4. * local.array() * (local.array() * alpha_squared_ + alpha_beta_) +
61 const Eigen::VectorXd jTj_inv = jTj.cwiseInverse().transpose();
62
63 return jacobian.array().rowwise() * jTj_inv.transpose().array();
64}
65
67 const Eigen::MatrixXd& local) const {
68 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
69 "local coordinates out of bounds for reference element");
70
71 const auto jTj =
72 4. * local.array() * (local.array() * alpha_squared_ + alpha_beta_) +
74
75 return jTj.cwiseSqrt().transpose();
76}
77
78std::unique_ptr<Geometry> SegmentO2::SubGeometry(dim_t codim, dim_t i) const {
79 if (codim == 0) {
80 LF_ASSERT_MSG(i == 0, "i is out of bounds");
81 return std::make_unique<SegmentO2>(coords_);
82 }
83 if (codim == 1) {
84 LF_ASSERT_MSG(0 <= i && i <= 2, "i is out of bounds");
85 return std::make_unique<Point>(coords_.col(i));
86 }
87 LF_VERIFY_MSG(false, "codim is out of bounds");
88}
89
90std::vector<std::unique_ptr<Geometry>> SegmentO2::ChildGeometry(
91 const RefinementPattern& ref_pat, base::dim_t codim) const {
92 LF_VERIFY_MSG(ref_pat.RefEl() == lf::base::RefEl::kSegment(),
93 "Refinement pattern for " << ref_pat.RefEl().ToString());
94 LF_VERIFY_MSG(codim < 2, "Illegal codim = " << codim);
95
96 const double hLattice = 1. / static_cast<double>(ref_pat.LatticeConst());
97 std::vector<std::unique_ptr<Geometry>> childGeoPtrs = {};
98
99 // get coordinates of childGeometries
100 std::vector<Eigen::Matrix<int, Eigen::Dynamic, Eigen::Dynamic>> childPolygons(
101 ref_pat.ChildPolygons(codim));
102
103 const size_t noChildren = childPolygons.size();
104 LF_VERIFY_MSG(
105 noChildren == ref_pat.NumChildren(codim),
106 "NumChildren " << noChildren << " <-> " << ref_pat.NumChildren(codim));
107
108 // create a geometry object for each child
109 for (size_t child = 0; child < noChildren; ++child) {
110 // codim == 0: a child segment is described by a polygon with three vertices
111 // codim == 1: a child point by a single point ("polygon with one corner")
112 LF_VERIFY_MSG(
113 childPolygons[child].rows() == 1,
114 "childPolygons[child].rows() = " << childPolygons[child].rows());
115 LF_VERIFY_MSG(
116 childPolygons[child].cols() == (2 - codim),
117 "childPolygons[child].cols() = " << childPolygons[child].cols());
118
119 Eigen::MatrixXd locCoords(1, 3 - 2 * codim);
120
121 switch (codim) {
122 case 0: {
123 // SegmentO2 requires the coordinate of the midpoint
124 locCoords << hLattice * childPolygons[child].cast<double>(),
125 (hLattice * childPolygons[child].cast<double>()).sum() / 2;
126 childGeoPtrs.push_back(std::make_unique<SegmentO2>(Global(locCoords)));
127
128 break;
129 }
130 case 1: {
131 locCoords << hLattice * childPolygons[child].cast<double>();
132 childGeoPtrs.push_back(std::make_unique<Point>(Global(locCoords)));
133
134 break;
135 }
136 default: {
137 LF_VERIFY_MSG(false, "Illegal co-dimension");
138 }
139 }
140 }
141
142 return childGeoPtrs;
143}
144
145} // namespace lf::geometry
static constexpr RefEl kSegment()
Returns the (1-dimensional) reference segment.
Definition ref_el.h:153
std::string ToString() const
Return a string representation of this Reference element.
Definition ref_el.h:458
base::RefEl::dim_t dim_t
Abstract interface class for encoding topological local refinement
lf::base::size_type LatticeConst() const
Provides information about lattice constant used.
virtual lf::base::size_type NumChildren(lf::base::dim_t codim) const =0
provide number of child entities of a given co-dimension to be created by refinement
virtual std::vector< Eigen::Matrix< int, Eigen::Dynamic, Eigen::Dynamic > > ChildPolygons(lf::base::dim_t codim) const =0
provide lattice reference coordinates of vertices of child polygons
lf::base::RefEl RefEl() const
Returns topological type of entity for which the current object is set up.
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
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