14#include "segment_o2.h"
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) {
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);
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);
51 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
52 "local coordinates out of bounds for reference element");
55 return ((
beta_ * local) + (
gamma_ * local.array().square().matrix()) +
56 (
delta_ * local.row(0).cwiseProduct(local.row(1))))
65 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
66 "local coordinates out of bounds for reference element");
68 Eigen::MatrixXd tmp(
gamma_.rows(), 2 * local.cols());
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();
75 Eigen::MatrixXd local_reversed = local.colwise().reverse();
76 local_reversed.resize(1, local.size());
78 return beta_.replicate(1, local.cols()) + tmp +
79 (local_reversed.replicate(
delta_.rows(), 1).array().colwise() *
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");
92 Eigen::MatrixXd jac =
Jacobian(local);
93 Eigen::MatrixXd jacInvGram(jac.rows(), jac.cols());
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();
103 for (
long i = 0; i < local.cols(); ++i) {
104 auto jacobian = jac.block(0, 2 * i, jac.rows(), 2);
106 auto A = jacobian.col(0);
107 auto B = jacobian.col(1);
110 Eigen::MatrixXd tmp(2, 2);
111 tmp << B.dot(B), -AB, -AB, A.dot(A);
113 jacInvGram.block(0, 2 * i, jac.rows(), 2) =
114 jacobian * tmp / tmp.determinant();
124 LF_VERIFY_MSG((0. <= local.array()).all() && (local.array() <= 1.).all(),
125 "local coordinates out of bounds for reference element");
127 Eigen::MatrixXd jac =
Jacobian(local);
128 Eigen::VectorXd intElem(local.cols());
131 for (
long i = 0; i < local.cols(); ++i) {
132 intElem(i) = std::abs(jac.block(0, 2 * i, 2, 2).determinant());
135 for (
long i = 0; i < local.cols(); ++i) {
136 auto jacobian = jac.block(0, 2 * i, jac.rows(), 2);
138 auto A = jacobian.col(0);
139 auto B = jacobian.col(1);
142 intElem(i) = std::sqrt(std::abs(A.dot(A) * B.dot(B) - AB * AB));
152 LF_ASSERT_MSG(i == 0,
"i is out of bounds");
153 return std::make_unique<TriaO2>(
coords_);
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)
164 LF_ASSERT_MSG(0 <= i && i <= 5,
"i is out of bounds");
165 return std::make_unique<Point>(
coords_.col(i));
168 LF_VERIFY_MSG(
false,
"codim is out of bounds")
177 LF_VERIFY_MSG(codim < 3,
"Illegal codim " << codim);
179 const double hLattice = 1. /
static_cast<double>(ref_pat.
LatticeConst());
180 std::vector<std::unique_ptr<Geometry>> childGeoPtrs = {};
183 std::vector<Eigen::Matrix<int, Eigen::Dynamic, Eigen::Dynamic>> childPolygons(
189 "NumChildren " << noChildren <<
" <-> " << ref_pat.
NumChildren(codim));
192 for (
size_t child = 0; child < noChildren; ++child) {
198 childPolygons[child].rows() == 2,
199 "childPolygons[child].rows() = " << childPolygons[child].rows());
201 childPolygons[child].cols() == (3 - codim),
202 "childPolygons[child].cols() = " << childPolygons[child].cols());
204 const Eigen::MatrixXd locChildPolygonCoords(
205 hLattice * childPolygons[child].cast<double>());
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.;
215 childGeoPtrs.push_back(
216 std::make_unique<TriaO2>(
Global(locChildCoords)));
221 Eigen::MatrixXd locChildCoords(locChildPolygonCoords.rows(), 3);
222 locChildCoords << locChildPolygonCoords,
223 (locChildPolygonCoords.col(0) + locChildPolygonCoords.col(1)) / 2.;
225 childGeoPtrs.push_back(
226 std::make_unique<SegmentO2>(
Global(locChildCoords)));
231 childGeoPtrs.push_back(
232 std::make_unique<Point>(
Global(locChildPolygonCoords)));
237 LF_VERIFY_MSG(
false,
"Illegal co-dimension");
static constexpr RefEl kTria()
Returns the reference triangle.
std::string ToString() const
Return a string representation of this Reference element.
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.
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Eigen::Matrix< double, Eigen::Dynamic, 6 > coords_
Coordinates of the 6 vertices/midpoints, stored in matrix columns.
Eigen::Matrix< double, Eigen::Dynamic, 2 > beta_
Eigen::MatrixXd JacobianInverseGramian(const Eigen::MatrixXd &local) const override
Evaluate the Jacobian * Inverse Gramian ( ) simultaneously at numPoints.
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::Matrix< double, Eigen::Dynamic, 1 > alpha_
Eigen::Matrix< double, Eigen::Dynamic, 2 > gamma_
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.
Eigen::Matrix< double, Eigen::Dynamic, 1 > delta_
Eigen::VectorXd IntegrationElement(const Eigen::MatrixXd &local) const override
The integration element (factor appearing in integral transformation formula, see below) at number of...
Eigen::Matrix< double, Eigen::Dynamic, 2 > gamma_x_2_
TriaO2(Eigen::Matrix< double, Eigen::Dynamic, 6 > coords)
Constructor building triangle from vertex/midpoint coordinates.
Eigen::MatrixXd Jacobian(const Eigen::MatrixXd &local) const override
Evaluate the jacobian of the mapping simultaneously at numPoints points.
unsigned int size_type
general type for variables related to size of arrays
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...