LehrFEM++ 1.0.0
A simple Finite Element Library for teaching
Loading...
Searching...
No Matches
tria_o1.cc
1#include "tria_o1.h"
2
3#include <Eigen/Eigen>
4
5#include "point.h"
6#include "segment_o1.h"
7
8namespace lf::geometry {
9
11 const Eigen::Matrix<double, Eigen::Dynamic, 3>& coords, double tol) {
12 // World dimension
13 const Geometry::dim_t wd = coords.rows();
14 // Length tests
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();
18 // Test lengths of edges versus circumference.
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");
23 // Area test
24 switch (wd) {
25 case 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);
31 return true;
32 break;
33 }
34 case 3: {
35 const Eigen::Matrix<double, 3, 3> c3d(coords.block<3, 3>(0, 0));
36 const double area =
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");
39 return true;
40 break;
41 }
42 default: {
43 LF_ASSERT_MSG(false, "Illegal world dimension" << wd);
44 break;
45 }
46 }
47 return false;
48}
49
50TriaO1::TriaO1(Eigen::Matrix<double, Eigen::Dynamic, 3> coords)
51 : coords_(std::move(coords)),
52 jacobian_(coords_.rows(), 2),
53 jacobian_inverse_gramian_(coords_.rows(), 2) {
54 // Make sure that the triangle has a proper shape.
56 // Precompute constant Jacobian
57 jacobian_ << coords_.col(1) - coords_.col(0), coords_.col(2) - coords_.col(0);
58 // Precompute constant metric factor
59 if (coords_.rows() == 2) {
60 jacobian_inverse_gramian_ = jacobian_.transpose().inverse();
61 integrationElement_ = std::abs(jacobian_.determinant());
62 } else {
63 jacobian_inverse_gramian_ = Eigen::MatrixXd(
64 jacobian_ * (jacobian_.transpose() * jacobian_).inverse());
66 std::sqrt((jacobian_.transpose() * jacobian_).determinant());
67 }
68} // end constructor
69
70Eigen::MatrixXd TriaO1::Global(const Eigen::MatrixXd& local) const {
71 return coords_.col(0) *
72 (1 - local.array().row(0) - local.array().row(1)).matrix() +
73 coords_.col(1) * local.row(0) + coords_.col(2) * local.row(1);
74}
75
76std::unique_ptr<Geometry> TriaO1::SubGeometry(dim_t codim, dim_t i) const {
77 using std::make_unique;
78 switch (codim) {
79 case 0:
80 LF_ASSERT_MSG(i == 0, "codim = 0: i = " << i << " is out of bounds.");
81 return std::make_unique<TriaO1>(coords_);
82 case 1:
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)),
88 coords_.col(RefEl().SubSubEntity2SubEntity(1, i, 1, 1)))
89 .finished());
90 case 2:
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));
94 default:
95 LF_VERIFY_MSG(false, "codim " << codim << " is out of bounds.");
96 }
97}
98
99std::vector<std::unique_ptr<Geometry>> TriaO1::ChildGeometry(
100 const RefinementPattern& ref_pat, base::dim_t codim) const {
101 // The refinement pattern must be for a triangle
102 LF_VERIFY_MSG(ref_pat.RefEl() == lf::base::RefEl::kTria(),
103 "Refinement pattern for " << ref_pat.RefEl().ToString());
104 LF_VERIFY_MSG(codim < 3, "Illegal codim " << codim);
105 // Lattice meshwidth
106 const double h_lattice = 1.0 / static_cast<double>(ref_pat.LatticeConst());
107 // Obtain geometry of children as lattice polygon
108 std::vector<Eigen::Matrix<int, Eigen::Dynamic, Eigen::Dynamic>>
109 child_polygons(ref_pat.ChildPolygons(codim));
110 // Number of child segments
111 const auto no_children = base::narrow<base::size_type>(child_polygons.size());
112 LF_VERIFY_MSG(
113 no_children == ref_pat.NumChildren(codim),
114 "no_children = " << no_children << " <-> " << ref_pat.NumChildren(codim));
115 // return variable
116 std::vector<std::unique_ptr<Geometry>> child_geo_uptrs{};
117 // For each child triangle create a geometry object and a unique pointer to
118 // it.
119 for (int l = 0; l < no_children; l++) {
120 // A single child triangle is described by a lattice polygon with
121 // three vertices, a child segment by a polygon with two vertices,
122 // a child point by a single point = "polygon with one corner"
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());
127 // Normalize lattice coordinates
128 const Eigen::MatrixXd child_geo(
129 Global(h_lattice * child_polygons[l].cast<double>()));
130 switch (codim) {
131 case 0: {
132 // child is a triangle
133 child_geo_uptrs.push_back(std::make_unique<TriaO1>(child_geo));
134 break;
135 }
136 case 1: {
137 // child is an edge
138 child_geo_uptrs.push_back(std::make_unique<SegmentO1>(child_geo));
139 break;
140 }
141 case 2: {
142 // child is a node
143 child_geo_uptrs.push_back(std::make_unique<Point>(child_geo));
144 break;
145 }
146 default:
147 LF_VERIFY_MSG(false, "codim out of bounds.");
148 } // end switch codim
149 } // end loop over the children
150 return (child_geo_uptrs);
151}
152
153} // 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.
TriaO1(Eigen::Matrix< double, Eigen::Dynamic, 3 > coords)
Definition tria_o1.cc:50
Eigen::MatrixXd Global(const Eigen::MatrixXd &local) const override
Map a number of points in local coordinates into the global coordinate system.
Definition tria_o1.cc:70
dim_t DimGlobal() const override
Dimension of the image of this mapping.
Definition tria_o1.h:28
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_
Definition tria_o1.h:63
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
Definition tria_o1.cc:99
Eigen::Matrix< double, Eigen::Dynamic, 3 > coords_
Definition tria_o1.h:62
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_o1.cc:76
double integrationElement_
Definition tria_o1.h:65
Eigen::Matrix< double, Eigen::Dynamic, 2 > jacobian_inverse_gramian_
Definition tria_o1.h:64
base::RefEl RefEl() const override
The Reference element that defines the domain of this mapping.
Definition tria_o1.h:29
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
bool assertNonDegenerateTriangle(const Eigen::Matrix< double, Eigen::Dynamic, 3 > &coords, double tol)
Asserting non-degenerate shape of a flat triangle.
Definition tria_o1.cc:10