LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
tria_o2.cc
1
9#include "tria_o2.h"
10
11#include <Eigen/Eigen>
12
13#include "point.h"
14#include "segment_o2.h"
15
16namespace lf::geometry {
17
18TriaO2::TriaO2(Eigen::Matrix<double, Eigen::Dynamic, 6> coords)
19 : coords_(std::move(coords)),
20 alpha_(coords_.rows()),
21 beta_(coords_.rows(), 2),
22 gamma_(coords_.rows(), 2),
23 delta_(coords_.rows()),
24 gamma_x_2_(coords_.rows(), 2) {
25 /*
26 * 2 C
27 * | \ | \
28 * 5 4 -> F E
29 * | \ | \
30 * 0 - 3 - 1 A - D - B
31 */
32 const Eigen::VectorXd& A = coords_.col(0);
33 const Eigen::VectorXd& B = coords_.col(1);
34 const Eigen::VectorXd& C = coords_.col(2);
35 const Eigen::VectorXd& D = coords_.col(3);
36 const Eigen::VectorXd& E = coords_.col(4);
37 const Eigen::VectorXd& F = coords_.col(5);
38
39 // Compute monomial coeffcients of componentwise quadratic
40 // polynomial
41 alpha_ << A;
42 beta_ << 4. * D - 3. * A - B, 4. * F - 3. * A - C;
43 gamma_ << 2. * (A + B) - 4. * D, 2. * (A + C) - 4. * F;
44 delta_ << 4. * (A + E - D - F);
45 // coefficient for Jacobian()
46 gamma_x_2_ << 2. * gamma_;
47}
48
49/* SAM_LISTING_BEGIN_1 */
50Eigen::MatrixXd TriaO2::Global(const Eigen::MatrixXd& local) const {
51 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
52 "local coordinates out of bounds for reference element");
53 // Direct vectorized evaluation of componentwise quadratic mapping
54 // given through its monomial coefficients.
55 return ((beta_ * local) + (gamma_ * local.array().square().matrix()) +
56 (delta_ * local.row(0).cwiseProduct(local.row(1))))
57 .colwise() +
58 alpha_;
59}
60
61/* SAM_LISTING_END_1 */
62
63/* SAM_LISTING_BEGIN_2 */
64Eigen::MatrixXd TriaO2::Jacobian(const Eigen::MatrixXd& local) const {
65 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
66 "local coordinates out of bounds for reference element");
67
68 Eigen::MatrixXd tmp(gamma_.rows(), 2 * local.cols());
69
70 for (long i = 0; i < local.cols(); ++i) {
71 tmp.block(0, 2 * i, tmp.rows(), 2) =
72 gamma_x_2_.array().rowwise() * local.col(i).transpose().array();
73 }
74
75 Eigen::MatrixXd local_reversed = local.colwise().reverse();
76 local_reversed.resize(1, local.size());
77
78 return beta_.replicate(1, local.cols()) + tmp +
79 (local_reversed.replicate(delta_.rows(), 1).array().colwise() *
80 delta_.array())
81 .matrix();
82}
83
84/* SAM_LISTING_END_2 */
85
86/* SAM_LISTING_BEGIN_3 */
88 const Eigen::MatrixXd& local) const {
89 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
90 "local coordinates out of bounds for reference element");
91
92 Eigen::MatrixXd jac = Jacobian(local);
93 Eigen::MatrixXd jacInvGram(jac.rows(), jac.cols());
94
95 if (DimGlobal() == 2) {
96 for (long i = 0; i < local.cols(); ++i) {
97 auto jacobian = jac.block(0, 2 * i, 2, 2);
98 jacInvGram.block(0, 2 * i, 2, 2) << jacobian(1, 1), -jacobian(1, 0),
99 -jacobian(0, 1), jacobian(0, 0);
100 jacInvGram.block(0, 2 * i, 2, 2) /= jacobian.determinant();
101 }
102 } else {
103 for (long i = 0; i < local.cols(); ++i) {
104 auto jacobian = jac.block(0, 2 * i, jac.rows(), 2);
105
106 auto A = jacobian.col(0);
107 auto B = jacobian.col(1);
108 auto AB = A.dot(B);
109
110 Eigen::MatrixXd tmp(2, 2);
111 tmp << B.dot(B), -AB, -AB, A.dot(A);
112
113 jacInvGram.block(0, 2 * i, jac.rows(), 2) =
114 jacobian * tmp / tmp.determinant();
115 }
116 }
117
118 return jacInvGram;
119}
120/* SAM_LISTING_END_3 */
121
122/* SAM_LISTING_BEGIN_4 */
123Eigen::VectorXd TriaO2::IntegrationElement(const Eigen::MatrixXd& local) const {
124 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
125 "local coordinates out of bounds for reference element");
126
127 Eigen::MatrixXd jac = Jacobian(local);
128 Eigen::VectorXd intElem(local.cols());
129
130 if (DimGlobal() == 2) {
131 for (long i = 0; i < local.cols(); ++i) {
132 intElem(i) = std::abs(jac.block(0, 2 * i, 2, 2).determinant());
133 }
134 } else {
135 for (long i = 0; i < local.cols(); ++i) {
136 auto jacobian = jac.block(0, 2 * i, jac.rows(), 2);
137
138 auto A = jacobian.col(0);
139 auto B = jacobian.col(1);
140 auto AB = A.dot(B);
141
142 intElem(i) = std::sqrt(std::abs(A.dot(A) * B.dot(B) - AB * AB));
143 }
144 }
145 return intElem;
146}
147/* SAM_LISTING_END_4 */
148
149std::unique_ptr<Geometry> TriaO2::SubGeometry(dim_t codim, dim_t i) const {
150 switch (codim) {
151 case 0: {
152 LF_ASSERT_MSG(i == 0, "i is out of bounds");
153 return std::make_unique<TriaO2>(coords_);
154 }
155 case 1: {
156 LF_ASSERT_MSG(0 <= i && i <= 2, "i is out of bounds");
157 return std::make_unique<SegmentO2>(
158 (Eigen::Matrix<double, Eigen::Dynamic, 3>(DimGlobal(), 3)
159 << coords_.col(i),
160 coords_.col((i + 1) % 3), coords_.col(i + 3))
161 .finished());
162 }
163 case 2: {
164 LF_ASSERT_MSG(0 <= i && i <= 5, "i is out of bounds");
165 return std::make_unique<Point>(coords_.col(i));
166 }
167 default: {
168 LF_VERIFY_MSG(false, "codim is out of bounds")
169 }
170 }
171}
172
173std::vector<std::unique_ptr<Geometry>> TriaO2::ChildGeometry(
174 const RefinementPattern& ref_pat, lf::base::dim_t codim) const {
175 LF_VERIFY_MSG(ref_pat.RefEl() == lf::base::RefEl::kTria(),
176 "Refinement pattern for " << ref_pat.RefEl().ToString());
177 LF_VERIFY_MSG(codim < 3, "Illegal codim " << codim);
178
179 const double hLattice = 1. / static_cast<double>(ref_pat.LatticeConst());
180 std::vector<std::unique_ptr<Geometry>> childGeoPtrs = {};
181
182 // get coordinates of childGeometries
183 std::vector<Eigen::Matrix<int, Eigen::Dynamic, Eigen::Dynamic>> childPolygons(
184 ref_pat.ChildPolygons(codim));
185
186 const base::size_type noChildren = childPolygons.size();
187 LF_VERIFY_MSG(
188 noChildren == ref_pat.NumChildren(codim),
189 "NumChildren " << noChildren << " <-> " << ref_pat.NumChildren(codim));
190
191 // create a geometry object for each child
192 for (size_t child = 0; child < noChildren; ++child) {
193 // codim == 0: a child triangle is described by a lattice polygon with six
194 // vertices
195 // codim == 1: a child segment is described by a polygon with three vertices
196 // codim == 2: a child point by a single point ("polygon with one corner")
197 LF_VERIFY_MSG(
198 childPolygons[child].rows() == 2,
199 "childPolygons[child].rows() = " << childPolygons[child].rows());
200 LF_VERIFY_MSG(
201 childPolygons[child].cols() == (3 - codim),
202 "childPolygons[child].cols() = " << childPolygons[child].cols());
203
204 const Eigen::MatrixXd locChildPolygonCoords(
205 hLattice * childPolygons[child].cast<double>());
206
207 switch (codim) {
208 case 0: {
209 Eigen::MatrixXd locChildCoords(locChildPolygonCoords.rows(), 6);
210 locChildCoords << locChildPolygonCoords,
211 (locChildPolygonCoords.col(0) + locChildPolygonCoords.col(1)) / 2.,
212 (locChildPolygonCoords.col(1) + locChildPolygonCoords.col(2)) / 2.,
213 (locChildPolygonCoords.col(2) + locChildPolygonCoords.col(0)) / 2.;
214
215 childGeoPtrs.push_back(
216 std::make_unique<TriaO2>(Global(locChildCoords)));
217
218 break;
219 }
220 case 1: {
221 Eigen::MatrixXd locChildCoords(locChildPolygonCoords.rows(), 3);
222 locChildCoords << locChildPolygonCoords,
223 (locChildPolygonCoords.col(0) + locChildPolygonCoords.col(1)) / 2.;
224
225 childGeoPtrs.push_back(
226 std::make_unique<SegmentO2>(Global(locChildCoords)));
227
228 break;
229 }
230 case 2: {
231 childGeoPtrs.push_back(
232 std::make_unique<Point>(Global(locChildPolygonCoords)));
233
234 break;
235 }
236 default: {
237 LF_VERIFY_MSG(false, "Illegal co-dimension");
238 }
239 }
240 }
241
242 return childGeoPtrs;
243}
244
245} // namespace lf::geometry
static constexpr RefEl kTria()
Returns the reference triangle.
Definition ref_el.h:161
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::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition tria_o2.cc:50
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition tria_o2.h:42
Eigen::Matrix< double, Eigen::Dynamic, 6 > coords_
Coordinates of the 6 vertices/midpoints, stored in matrix columns.
Definition tria_o2.h:66
Eigen::Matrix< double, Eigen::Dynamic, 2 > beta_
Definition tria_o2.h:73
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition tria_o2.cc:87
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_o2.cc:149
Eigen::Matrix< double, Eigen::Dynamic, 1 > alpha_
Definition tria_o2.h:72
Eigen::Matrix< double, Eigen::Dynamic, 2 > gamma_
Definition tria_o2.h:74
std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pat, lf::base::dim_t codim) const override
Generate geometry objects for child entities created in the course of refinement.
Definition tria_o2.cc:173
Eigen::Matrix< double, Eigen::Dynamic, 1 > delta_
Definition tria_o2.h:75
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_o2.cc:123
Eigen::Matrix< double, Eigen::Dynamic, 2 > gamma_x_2_
Definition tria_o2.h:78
TriaO2(Eigen::Matrix< double, Eigen::Dynamic, 6 > coords)
Constructor building triangle from vertex/midpoint coordinates.
Definition tria_o2.cc:18
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition tria_o2.cc:64
unsigned int size_type
general type for variables related to size of arrays
Definition types.h:20
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