14#ifndef DUMUX_GEOMETRY_TRIANGULATION_HH
15#define DUMUX_GEOMETRY_TRIANGULATION_HH
23#include <dune/common/exceptions.hh>
24#include <dune/common/fvector.hh>
31namespace TriangulationPolicy {
46using DefaultDimPolicies = std::tuple<DelaunayPolicy, MidPointPolicy, ConvexHullPolicy>;
51template<
int dim,
int dimWorld>
52using DefaultPolicy = std::tuple_element_t<dim-1, Detail::DefaultDimPolicies>;
63template<
int dim,
int dimWorld,
class ctype>
64using Triangulation = std::vector< std::array<Dune::FieldVector<ctype, dimWorld>, dim+1> >;
76template<
int dim,
int dimWorld,
class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>,
77 class RandomAccessContainer,
78 std::enable_if_t< std::is_same_v<Policy, TriangulationPolicy::DelaunayPolicy>
79 && dim == 1,
int> = 0 >
83 using ctype =
typename RandomAccessContainer::value_type::value_type;
84 using Point = Dune::FieldVector<ctype, dimWorld>;
86 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
87 "Triangulation expects Dune::FieldVector as point type");
89 if (points.size() == 2)
93 assert(points.size() > 1);
94 DUNE_THROW(Dune::NotImplemented,
"1d triangulation for point cloud size > 2");
109template<
int dim,
int dimWorld,
class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>,
110 class RandomAccessContainer,
111 std::enable_if_t< std::is_same_v<Policy, TriangulationPolicy::M
idPo
intPolicy>
112 && dim == 2,
int> = 0 >
113inline Triangulation<dim, dimWorld, typename RandomAccessContainer::value_type::value_type>
114triangulate(
const RandomAccessContainer& convexHullPoints)
116 using ctype =
typename RandomAccessContainer::value_type::value_type;
117 using Point = Dune::FieldVector<ctype, dimWorld>;
118 using Triangle = std::array<Point, 3>;
120 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
121 "Triangulation expects Dune::FieldVector as point type");
123 if (convexHullPoints.size() < 3)
124 DUNE_THROW(Dune::InvalidStateException,
"Try to triangulate point cloud with less than 3 points!");
126 if (convexHullPoints.size() == 3)
127 return std::vector<Triangle>(1, {convexHullPoints[0], convexHullPoints[1], convexHullPoints[2]});
130 for (
const auto& p : convexHullPoints)
132 midPoint /= convexHullPoints.size();
134 std::vector<Triangle> triangulation;
135 triangulation.reserve(convexHullPoints.size());
137 for (std::size_t i = 0; i < convexHullPoints.size()-1; ++i)
138 triangulation.emplace_back(Triangle{midPoint, convexHullPoints[i], convexHullPoints[i+1]});
140 triangulation.emplace_back(Triangle{midPoint, convexHullPoints[convexHullPoints.size()-1], convexHullPoints[0]});
142 return triangulation;
157template<
int dim,
int dimWorld,
class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>,
158 class RandomAccessContainer,
159 std::enable_if_t< std::is_same_v<Policy, TriangulationPolicy::ConvexHullPolicy>
160 && dim == 2,
int> = 0 >
161inline Triangulation<dim, dimWorld, typename RandomAccessContainer::value_type::value_type>
164 const auto convexHullPoints = grahamConvexHull<2>(points);
165 return triangulate<dim, dimWorld, TriangulationPolicy::MidPointPolicy>(convexHullPoints);
178template<
int dim,
int dimWorld,
class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>,
179 class RandomAccessContainer,
180 std::enable_if_t< std::is_same_v<Policy, TriangulationPolicy::ConvexHullPolicy>
181 && dim == 3,
int> = 0 >
182inline Triangulation<dim, dimWorld, typename RandomAccessContainer::value_type::value_type>
185 using ctype =
typename RandomAccessContainer::value_type::value_type;
186 using Point = Dune::FieldVector<ctype, dimWorld>;
187 using Tetrahedron = std::array<Point, 4>;
189 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
190 "Triangulation expects Dune::FieldVector as point type");
192 const auto numPoints = points.size();
194 DUNE_THROW(Dune::InvalidStateException,
"Trying to create 3D triangulation of point cloud with less than 4 points!");
197 return std::vector<Tetrahedron>(1, {points[0], points[1], points[2], points[3]});
202 Point lowerLeft(1e100);
203 Point upperRight(-1e100);
204 for (
const auto& p : points)
207 for (
int i = 0; i < dimWorld; ++i)
209 using std::max;
using std::min;
210 lowerLeft[i] = min(p[i], lowerLeft[i]);
211 upperRight[i] = max(p[i], upperRight[i]);
214 midPoint /= numPoints;
216 auto magnitude = 0.0;
218 for (
int i = 0; i < dimWorld; ++i)
219 magnitude = max(upperRight[i] - lowerLeft[i], magnitude);
220 const auto eps = 1e-7*magnitude;
221 const auto eps2 = eps*magnitude;
224 std::vector<Tetrahedron> triangulation;
225 triangulation.reserve(numPoints);
228 std::vector<Point> coplanarPointBuffer;
229 coplanarPointBuffer.reserve(std::min<std::size_t>(12, numPoints-1));
234 std::vector<std::pair<int, Point>> coplanarClusters;
235 coplanarClusters.reserve(numPoints/3);
241 for (
int i = 0; i < numPoints; ++i)
243 for (
int j = i+1; j < numPoints; ++j)
245 for (
int k = j+1; k < numPoints; ++k)
247 const auto pointI = points[i];
248 const auto ab = points[j] - pointI;
249 const auto ac = points[k] - pointI;
253 coplanarPointBuffer.clear();
257 const bool isAdmissible = [&]()
261 if (
normal.two_norm2() < eps2*eps2)
265 for (
int m = 0; m < numPoints; ++m)
267 if (m != i && m != j && m != k)
270 const auto ad = points[m] - pointI;
271 const auto sp =
normal*ad;
274 using std::abs;
using std::signbit;
275 const bool coplanar = abs(sp) < eps2*magnitude;
276 int newMarker = coplanar ? 0 : signbit(sp) ? -1 : 1;
280 if (marker == 0 && newMarker != 0)
285 if (newMarker != 0 && marker != newMarker)
292 if (m < k && std::find_if(
293 coplanarClusters.begin(), coplanarClusters.end(),
294 [=](
const auto& c){ return c.first == std::min(m, i) && abs(ab*c.second) < eps2*magnitude; }
295 ) != coplanarClusters.end())
298 coplanarPointBuffer.clear();
302 coplanarPointBuffer.push_back(points[m]);
309 if (!coplanarPointBuffer.empty())
311 coplanarPointBuffer.insert(coplanarPointBuffer.end(), { points[i], points[j], points[k] });
312 coplanarClusters.emplace_back(std::make_pair(i,
normal));
326 if (!coplanarPointBuffer.empty())
328 const auto triangles = triangulate<2, 3, TriangulationPolicy::ConvexHullPolicy>(coplanarPointBuffer);
329 for (
const auto& triangle : triangles)
331 const auto ab = triangle[1] - triangle[0];
332 const auto ac = triangle[2] - triangle[0];
334 const auto am = midPoint - triangle[0];
335 const auto sp =
normal*am;
337 const bool isBelow = signbit(sp);
339 triangulation.emplace_back(Tetrahedron{
340 triangle[0], triangle[2], triangle[1], midPoint
343 triangulation.emplace_back(Tetrahedron{
344 triangle[0], triangle[1], triangle[2], midPoint
350 const auto am = midPoint - pointI;
351 const auto sp =
normal*am;
353 const bool isBelow = signbit(sp);
355 triangulation.emplace_back(Tetrahedron{
356 pointI, points[k], points[j], midPoint
359 triangulation.emplace_back(Tetrahedron{
360 pointI, points[j], points[k], midPoint
369 if (triangulation.size() < 4)
370 DUNE_THROW(Dune::InvalidStateException,
"Something went wrong with the triangulation!");
372 return triangulation;
A function to compute the convex hull of a point cloud.
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:642
Triangulation< dim, dimWorld, typename RandomAccessContainer::value_type::value_type > triangulate(const RandomAccessContainer &points)
Triangulate area given points of a convex hull (1d)
Definition: triangulation.hh:81
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
std::vector< std::array< Dune::FieldVector< ctype, dimWorld >, dim+1 > > Triangulation
The default data type to store triangulations.
Definition: triangulation.hh:64
Detect if a point intersects a simplex (including boundary)
Define some often used mathematical functions.
std::tuple_element_t< dim-1, Detail::DefaultDimPolicies > DefaultPolicy
Default policy for a given dimension.
Definition: triangulation.hh:52
Definition: triangulation.hh:39
Delaunay-type triangulations.
Definition: triangulation.hh:42
Definition: triangulation.hh:35