11 const Eigen::Matrix<double, Eigen::Dynamic, 3>& coords,
double tol) {
15 const double e0lensq = (coords.col(1) - coords.col(0)).squaredNorm();
16 const double e1lensq = (coords.col(2) - coords.col(1)).squaredNorm();
17 const double e2lensq = (coords.col(0) - coords.col(2)).squaredNorm();
19 const double circum = e0lensq + e1lensq + e2lensq;
20 LF_VERIFY_MSG(e0lensq > tol * circum,
"Collapsed edge 0");
21 LF_VERIFY_MSG(e1lensq > tol * circum,
"Collapsed edge 1");
22 LF_VERIFY_MSG(e2lensq > tol * circum,
"Collapsed edge 2");
26 const double area = std::fabs(
27 ((coords(0, 1) - coords(0, 0)) * (coords(1, 2) - coords(1, 0)) -
28 (coords(1, 1) - coords(1, 0)) * (coords(0, 2) - coords(0, 0))));
29 LF_VERIFY_MSG(area > tol * circum,
30 "Degenerate 2D triangle, geo = " << coords);
35 const Eigen::Matrix<double, 3, 3> c3d(coords.block<3, 3>(0, 0));
37 ((c3d.col(1) - c3d.col(0)).cross(c3d.col(2) - c3d.col(0))).norm();
38 LF_VERIFY_MSG(area > tol * circum,
"Degenerate 3D triangle");
43 LF_ASSERT_MSG(
false,
"Illegal world dimension" << wd);
51 : coords_(std::move(coords)),
52 jacobian_(coords_.rows(), 2),
53 jacobian_inverse_gramian_(coords_.rows(), 2) {
72 (1 - local.array().row(0) - local.array().row(1)).matrix() +
77 using std::make_unique;
80 LF_ASSERT_MSG(i == 0,
"codim = 0: i = " << i <<
" is out of bounds.");
81 return std::make_unique<TriaO1>(
coords_);
83 LF_ASSERT_MSG(i >= 0 && i < 3,
84 "codim =1: i = " << i <<
" is out of bounds.");
85 return make_unique<SegmentO1>(
86 (Eigen::Matrix<double, Eigen::Dynamic, 2>(
DimGlobal(), 2)
87 <<
coords_.col(
RefEl().SubSubEntity2SubEntity(1, i, 1, 0)),
91 LF_ASSERT_MSG(i >= 0 && i < 3,
92 "codim = 2: i = " << i <<
" is out of bounds.");
93 return make_unique<Point>(
coords_.col(i));
95 LF_VERIFY_MSG(
false,
"codim " << codim <<
" is out of bounds.");
104 LF_VERIFY_MSG(codim < 3,
"Illegal codim " << codim);
106 const double h_lattice = 1.0 /
static_cast<double>(ref_pat.
LatticeConst());
108 std::vector<Eigen::Matrix<int, Eigen::Dynamic, Eigen::Dynamic>>
111 const auto no_children = base::narrow<base::size_type>(child_polygons.size());
114 "no_children = " << no_children <<
" <-> " << ref_pat.
NumChildren(codim));
116 std::vector<std::unique_ptr<Geometry>> child_geo_uptrs{};
119 for (
int l = 0; l < no_children; l++) {
123 LF_VERIFY_MSG(child_polygons[l].rows() == 2,
124 "child_polygons[l].rows() = " << child_polygons[l].rows());
125 LF_VERIFY_MSG(child_polygons[l].cols() == 3 - codim,
126 "child_polygons[l].cols() = " << child_polygons[l].cols());
128 const Eigen::MatrixXd child_geo(
129 Global(h_lattice * child_polygons[l].cast<double>()));
133 child_geo_uptrs.push_back(std::make_unique<TriaO1>(child_geo));
138 child_geo_uptrs.push_back(std::make_unique<SegmentO1>(child_geo));
143 child_geo_uptrs.push_back(std::make_unique<Point>(child_geo));
147 LF_VERIFY_MSG(
false,
"codim out of bounds.");
150 return (child_geo_uptrs);
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.
TriaO1(Eigen::Matrix< double, Eigen::Dynamic, 3 > coords)
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, 2 > jacobian_
std::vector< std::unique_ptr< Geometry > > ChildGeometry(const RefinementPattern &ref_pat, lf::base::dim_t codim) const override
creation of child geometries as specified in refinement pattern
Eigen::Matrix< double, Eigen::Dynamic, 3 > coords_
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...
double integrationElement_
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_inverse_gramian_
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
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...
bool assertNonDegenerateTriangle(const Eigen::Matrix< double, Eigen::Dynamic, 3 > &coords, double tol)
Asserting non-degenerate shape of a flat triangle.