LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
quad_o2.cc
1
9#include "quad_o2.h"
10
11#include <Eigen/Eigen>
12
13#include "point.h"
14
15namespace lf::geometry {
16
17QuadO2::QuadO2(Eigen::Matrix<double, Eigen::Dynamic, 8> coords)
18 : coords_(std::move(coords)),
19 alpha_(coords_.rows()),
20 beta_(coords_.rows(), 2),
21 gamma_(coords_.rows(), 2),
22 delta_(coords_.rows()),
23 epsilon_(coords_.rows(), 2),
24 gamma_x_2_(coords_.rows(), 2),
25 epsilon_x_2_(coords_.rows(), 2) {
26 /*
27 * 3 - 6 - 2 D - G - C
28 * | | | |
29 * 7 5 -> H F
30 * | | | |
31 * 0 - 4 - 1 A - E - B
32 */
33 const Eigen::VectorXd& A = coords_.col(0);
34 const Eigen::VectorXd& B = coords_.col(1);
35 const Eigen::VectorXd& C = coords_.col(2);
36 const Eigen::VectorXd& D = coords_.col(3);
37 const Eigen::VectorXd& E = coords_.col(4);
38 const Eigen::VectorXd& F = coords_.col(5);
39 const Eigen::VectorXd& G = coords_.col(6);
40 const Eigen::VectorXd& H = coords_.col(7);
41
42 alpha_ << A;
43 beta_ << -3. * A - B + 4. * E, -3. * A - D + 4. * H;
44 gamma_ << 2. * (A + B) - 4. * E, 2. * (A + D) - 4. * H;
45 delta_ << 5. * A - B - 3. * C - D - 4. * (E - F - G + H);
46 epsilon_ << -2. * (A + B - C - D) + 4. * (E - G),
47 -2. * (A - B - C + D) - 4. * (F - H);
48
49 // coefficients for Jacobian()
50 gamma_x_2_ << 2. * gamma_;
51 epsilon_x_2_ << 2. * epsilon_;
52}
53
54Eigen::MatrixXd QuadO2::Global(const Eigen::MatrixXd& local) const {
55 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
56 "local coordinates out of bounds for reference element");
57
58 auto local_squared = local.array().square();
59 return ((beta_ * local) + (gamma_ * local_squared.matrix()) +
60 (delta_ * local.row(0).cwiseProduct(local.row(1))) +
61 (epsilon_ *
62 (local_squared * local.colwise().reverse().array()).matrix()))
63 .colwise() +
64 alpha_;
65}
66
67Eigen::MatrixXd QuadO2::Jacobian(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 auto local_squared_reversed = local.array().square().colwise().reverse();
72 auto local_colwise_product = local.row(0).cwiseProduct(local.row(1));
73
74 Eigen::MatrixXd tmp_epsilon(gamma_.rows(), 2 * local.cols());
75 Eigen::MatrixXd tmp_epsilon_x_2 = epsilon_x_2_.replicate(1, local.cols());
76 Eigen::MatrixXd tmp_gamma(gamma_.rows(), 2 * local.cols());
77
78 for (long i = 0; i < local.cols(); ++i) {
79 tmp_epsilon.block(0, 2 * i, tmp_epsilon.rows(), 2) =
80 epsilon_.rowwise().reverse().array().rowwise() *
81 local_squared_reversed.col(i).transpose().array();
82 tmp_epsilon_x_2.block(0, 2 * i, tmp_epsilon_x_2.rows(), 2) *=
83 local_colwise_product(i);
84 tmp_gamma.block(0, 2 * i, tmp_gamma.rows(), 2) =
85 gamma_x_2_.array().rowwise() * local.col(i).transpose().array();
86 }
87
88 Eigen::MatrixXd local_reversed = local.colwise().reverse();
89 local_reversed.resize(1, local.size());
90
91 return beta_.replicate(1, local.cols()) + tmp_gamma + tmp_epsilon +
92 tmp_epsilon_x_2 +
93 (local_reversed.replicate(delta_.rows(), 1).array().colwise() *
94 delta_.array())
95 .matrix();
96}
97
99 const Eigen::MatrixXd& local) const {
100 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
101 "local coordinates out of bounds for reference element");
102
103 Eigen::MatrixXd jac = Jacobian(local);
104 Eigen::MatrixXd jacInvGram(jac.rows(), jac.cols());
105
106 if (DimGlobal() == 2) {
107 for (long i = 0; i < local.cols(); ++i) {
108 auto jacobian = jac.block(0, 2 * i, 2, 2);
109 jacInvGram.block(0, 2 * i, 2, 2) << jacobian(1, 1), -jacobian(1, 0),
110 -jacobian(0, 1), jacobian(0, 0);
111 jacInvGram.block(0, 2 * i, 2, 2) /= jacobian.determinant();
112 }
113 } else {
114 for (long i = 0; i < local.cols(); ++i) {
115 auto jacobian = jac.block(0, 2 * i, jac.rows(), 2);
116
117 auto A = jacobian.col(0);
118 auto B = jacobian.col(1);
119 auto AB = A.dot(B);
120
121 Eigen::MatrixXd tmp(2, 2);
122 tmp << B.dot(B), -AB, -AB, A.dot(A);
123
124 jacInvGram.block(0, 2 * i, jac.rows(), 2) =
125 jacobian * tmp / tmp.determinant();
126 }
127 }
128
129 return jacInvGram;
130}
131
132Eigen::VectorXd QuadO2::IntegrationElement(const Eigen::MatrixXd& local) const {
133 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
134 "local coordinates out of bounds for reference element");
135
136 Eigen::MatrixXd jac = Jacobian(local);
137 Eigen::VectorXd intElem(local.cols());
138
139 if (DimGlobal() == 2) {
140 for (long i = 0; i < local.cols(); ++i) {
141 intElem(i) = std::abs(jac.block(0, 2 * i, 2, 2).determinant());
142 }
143 } else {
144 for (long i = 0; i < local.cols(); ++i) {
145 auto jacobian = jac.block(0, 2 * i, jac.rows(), 2);
146
147 auto A = jacobian.col(0);
148 auto B = jacobian.col(1);
149 auto AB = A.dot(B);
150
151 intElem(i) = std::sqrt(std::abs(A.dot(A) * B.dot(B) - AB * AB));
152 }
153 }
154
155 return intElem;
156}
157
158std::unique_ptr<Geometry> QuadO2::SubGeometry(dim_t codim, dim_t i) const {
159 switch (codim) {
160 case 0: {
161 LF_ASSERT_MSG(i == 0, "i is out of bounds");
162 return std::make_unique<QuadO2>(coords_);
163 }
164 case 1: {
165 LF_ASSERT_MSG(0 <= i && i <= 3, "i is out of bounds");
166 return std::make_unique<SegmentO2>(
167 (Eigen::Matrix<double, Eigen::Dynamic, 3>(DimGlobal(), 3)
168 << coords_.col(i),
169 coords_.col((i + 1) % 4), coords_.col(i + 4))
170 .finished());
171 }
172 case 2: {
173 LF_ASSERT_MSG(0 <= i && i <= 8, "i is out of bounds");
174 return std::make_unique<Point>(coords_.col(i));
175 }
176 default: {
177 LF_VERIFY_MSG(false, "codim is out of bounds")
178 }
179 }
180}
181
182std::vector<std::unique_ptr<Geometry>> QuadO2::ChildGeometry(
183 const RefinementPattern& ref_pat, lf::base::dim_t codim) const {
184 LF_VERIFY_MSG(ref_pat.RefEl() == lf::base::RefEl::kQuad(),
185 "Refinement pattern for " << ref_pat.RefEl().ToString());
186 LF_VERIFY_MSG(codim < 3, "Illegal codim " << codim);
187
188 const double hLattice = 1. / static_cast<double>(ref_pat.LatticeConst());
189 std::vector<std::unique_ptr<Geometry>> childGeoPtrs = {};
190
191 // get coordinates of childGeometries
192 std::vector<Eigen::Matrix<int, Eigen::Dynamic, Eigen::Dynamic>> childPolygons(
193 ref_pat.ChildPolygons(codim));
194
195 const int noChildren = childPolygons.size();
196 LF_VERIFY_MSG(
197 noChildren == ref_pat.NumChildren(codim),
198 "NumChildren " << noChildren << " <-> " << ref_pat.NumChildren(codim));
199
200 // create a geometry object for each child
201 for (size_t child = 0; child < noChildren; ++child) {
202 // codim == 0: a child triangle/quadrilateral is described by a lattice
203 // polygon with six/eight vertices
204 // codim == 1: a child segment is described by a polygon with three vertices
205 // codim == 2: a child point by a single point ("polygon with one corner")
206 LF_VERIFY_MSG(
207 childPolygons[child].rows() == 2,
208 "childPolygons[child].rows() = " << childPolygons[child].rows());
209 LF_VERIFY_MSG(
210 (childPolygons[child].cols() == (3 - codim)) ||
211 (childPolygons[child].cols() == 4),
212 "childPolygons[child].cols() = " << childPolygons[child].cols());
213
214 const Eigen::MatrixXd locChildPolygonCoords(
215 hLattice * childPolygons[child].cast<double>());
216
217 switch (codim) {
218 case 0: {
219 if (childPolygons[child].cols() == 3) {
220 Eigen::MatrixXd locChildCoords(locChildPolygonCoords.rows(), 6);
221 locChildCoords << locChildPolygonCoords,
222 (locChildPolygonCoords.col(0) + locChildPolygonCoords.col(1)) /
223 2.,
224 (locChildPolygonCoords.col(1) + locChildPolygonCoords.col(2)) /
225 2.,
226 (locChildPolygonCoords.col(2) + locChildPolygonCoords.col(0)) /
227 2.;
228
229 childGeoPtrs.push_back(
230 std::make_unique<TriaO2>(Global(locChildCoords)));
231
232 } else if (childPolygons[child].cols() == 4) {
233 Eigen::MatrixXd locChildCoords(locChildPolygonCoords.rows(), 8);
234 locChildCoords << locChildPolygonCoords,
235 (locChildPolygonCoords.col(0) + locChildPolygonCoords.col(1)) /
236 2.,
237 (locChildPolygonCoords.col(1) + locChildPolygonCoords.col(2)) /
238 2.,
239 (locChildPolygonCoords.col(2) + locChildPolygonCoords.col(3)) /
240 2.,
241 (locChildPolygonCoords.col(3) + locChildPolygonCoords.col(0)) /
242 2.;
243
244 childGeoPtrs.push_back(
245 std::make_unique<QuadO2>(Global(locChildCoords)));
246
247 } else {
248 LF_VERIFY_MSG(false, "childPolygons[" << child << "].cols() = "
249 << childPolygons[child].cols());
250 }
251
252 break;
253 }
254 case 1: {
255 Eigen::MatrixXd locChildCoords(locChildPolygonCoords.rows(), 3);
256 locChildCoords << locChildPolygonCoords,
257 (locChildPolygonCoords.col(0) + locChildPolygonCoords.col(1)) / 2.;
258
259 childGeoPtrs.push_back(
260 std::make_unique<SegmentO2>(Global(locChildCoords)));
261
262 break;
263 }
264 case 2: {
265 childGeoPtrs.push_back(
266 std::make_unique<Point>(Global(locChildPolygonCoords)));
267
268 break;
269 }
270 default: {
271 LF_VERIFY_MSG(false, "Illegal co-dimension");
272 }
273 }
274 }
275
276 return childGeoPtrs;
277}
278
279} // namespace lf::geometry
static constexpr RefEl kQuad()
Returns the reference quadrilateral.
Definition: ref_el.h:166
std::string ToString() const
Return a string representation of this Reference element.
Definition: ref_el.h:455
base::RefEl::dim_t dim_t
Eigen::Matrix< double, Eigen::Dynamic, 2 > gamma_x_2_
Definition: quad_o2.h:84
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
Definition: quad_o2.cc:67
Eigen::Matrix< double, Eigen::Dynamic, 2 > epsilon_x_2_
Definition: quad_o2.h:85
Eigen::Matrix< double, Eigen::Dynamic, 2 > beta_
Definition: quad_o2.h:78
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Definition: quad_o2.cc:132
Eigen::Matrix< double, Eigen::Dynamic, 2 > gamma_
Definition: quad_o2.h:79
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition: quad_o2.cc:54
Eigen::Matrix< double, Eigen::Dynamic, 8 > coords_
Coordinates of the 8 vertices/midpoints, stored in matrix columns.
Definition: quad_o2.h:70
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
Definition: quad_o2.cc:98
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: quad_o2.cc:182
Eigen::Matrix< double, Eigen::Dynamic, 1 > alpha_
Definition: quad_o2.h:77
Eigen::Matrix< double, Eigen::Dynamic, 1 > delta_
Definition: quad_o2.h:80
Eigen::Matrix< double, Eigen::Dynamic, 2 > epsilon_
Definition: quad_o2.h:81
QuadO2(Eigen::Matrix< double, Eigen::Dynamic, 8 > coords)
Constructor building quadrilateral from vertex/midpoint coordinates.
Definition: quad_o2.cc:17
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition: quad_o2.h:39
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: quad_o2.cc:158
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.
unsigned int dim_t
type for dimensions and co-dimensions and numbers derived from them
Definition: base.h:36
Defines the Geometry interface and provides a number of classes that implement this interface + addit...
Definition: compose.h:13