version 3.8
triangulation.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
14#ifndef DUMUX_GEOMETRY_TRIANGULATION_HH
15#define DUMUX_GEOMETRY_TRIANGULATION_HH
16
17#include <vector>
18#include <array>
19#include <algorithm>
20#include <type_traits>
21#include <tuple>
22
23#include <dune/common/exceptions.hh>
24#include <dune/common/fvector.hh>
25
26#include <dumux/common/math.hh>
29
30namespace Dumux {
31namespace TriangulationPolicy {
32
36
40
43
44#ifndef DOXYGEN
45namespace Detail {
46using DefaultDimPolicies = std::tuple<DelaunayPolicy, MidPointPolicy, ConvexHullPolicy>;
47} // end namespace Detail
48#endif
49
51template<int dim, int dimWorld>
52using DefaultPolicy = std::tuple_element_t<dim-1, Detail::DefaultDimPolicies>;
53
54} // end namespace TriangulationPolicy
55
63template<int dim, int dimWorld, class ctype>
64using Triangulation = std::vector< std::array<Dune::FieldVector<ctype, dimWorld>, dim+1> >;
65
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 >
81triangulate(const RandomAccessContainer& points)
82{
83 using ctype = typename RandomAccessContainer::value_type::value_type;
84 using Point = Dune::FieldVector<ctype, dimWorld>;
85
86 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
87 "Triangulation expects Dune::FieldVector as point type");
88
89 if (points.size() == 2)
90 return Triangulation<dim, dimWorld, ctype>({ {points[0], points[1]} });
91
93 assert(points.size() > 1);
94 DUNE_THROW(Dune::NotImplemented, "1d triangulation for point cloud size > 2");
95}
96
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::MidPointPolicy>
112 && dim == 2, int> = 0 >
113inline Triangulation<dim, dimWorld, typename RandomAccessContainer::value_type::value_type>
114triangulate(const RandomAccessContainer& convexHullPoints)
115{
116 using ctype = typename RandomAccessContainer::value_type::value_type;
117 using Point = Dune::FieldVector<ctype, dimWorld>;
118 using Triangle = std::array<Point, 3>;
119
120 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
121 "Triangulation expects Dune::FieldVector as point type");
122
123 if (convexHullPoints.size() < 3)
124 DUNE_THROW(Dune::InvalidStateException, "Try to triangulate point cloud with less than 3 points!");
125
126 if (convexHullPoints.size() == 3)
127 return std::vector<Triangle>(1, {convexHullPoints[0], convexHullPoints[1], convexHullPoints[2]});
128
129 Point midPoint(0.0);
130 for (const auto& p : convexHullPoints)
131 midPoint += p;
132 midPoint /= convexHullPoints.size();
133
134 std::vector<Triangle> triangulation;
135 triangulation.reserve(convexHullPoints.size());
136
137 for (std::size_t i = 0; i < convexHullPoints.size()-1; ++i)
138 triangulation.emplace_back(Triangle{midPoint, convexHullPoints[i], convexHullPoints[i+1]});
139
140 triangulation.emplace_back(Triangle{midPoint, convexHullPoints[convexHullPoints.size()-1], convexHullPoints[0]});
141
142 return triangulation;
143}
144
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>
162triangulate(const RandomAccessContainer& points)
163{
164 const auto convexHullPoints = grahamConvexHull<2>(points);
165 return triangulate<dim, dimWorld, TriangulationPolicy::MidPointPolicy>(convexHullPoints);
166}
167
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>
183triangulate(const RandomAccessContainer& points)
184{
185 using ctype = typename RandomAccessContainer::value_type::value_type;
186 using Point = Dune::FieldVector<ctype, dimWorld>;
187 using Tetrahedron = std::array<Point, 4>;
188
189 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
190 "Triangulation expects Dune::FieldVector as point type");
191
192 const auto numPoints = points.size();
193 if (numPoints < 4)
194 DUNE_THROW(Dune::InvalidStateException, "Trying to create 3D triangulation of point cloud with less than 4 points!");
195
196 if (numPoints == 4)
197 return std::vector<Tetrahedron>(1, {points[0], points[1], points[2], points[3]});
198
199 // compute the mid point of the point cloud (not the midpoint of the convex hull but this
200 // should not matter too much for the applications we have in mind here)
201 Point midPoint(0.0);
202 Point lowerLeft(1e100);
203 Point upperRight(-1e100);
204 for (const auto& p : points)
205 {
206 midPoint += p;
207 for (int i = 0; i < dimWorld; ++i)
208 {
209 using std::max; using std::min;
210 lowerLeft[i] = min(p[i], lowerLeft[i]);
211 upperRight[i] = max(p[i], upperRight[i]);
212 }
213 }
214 midPoint /= numPoints;
215
216 auto magnitude = 0.0;
217 using std::max;
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;
222
223 // reserve memory conservatively to avoid reallocation
224 std::vector<Tetrahedron> triangulation;
225 triangulation.reserve(numPoints);
226
227 // make a buffer for storing coplanar points and indices
228 std::vector<Point> coplanarPointBuffer;
229 coplanarPointBuffer.reserve(std::min<std::size_t>(12, numPoints-1));
230
231 // remember coplanar cluster planes
232 // coplanar clusters are uniquely identified by a point and a normal (plane)
233 // we only want to add each cluster once (when handling the first triangle in the cluster)
234 std::vector<std::pair<int, Point>> coplanarClusters;
235 coplanarClusters.reserve(numPoints/3);
236
237 // brute force algorithm: Try all possible triangles and check
238 // if they are triangles of the convex hull. This is achieved by checking if all
239 // other points are in the same half-space (side)
240 // the algorithm is O(n^4), so only do this for very small point clouds
241 for (int i = 0; i < numPoints; ++i)
242 {
243 for (int j = i+1; j < numPoints; ++j)
244 {
245 for (int k = j+1; k < numPoints; ++k)
246 {
247 const auto pointI = points[i];
248 const auto ab = points[j] - pointI;
249 const auto ac = points[k] - pointI;
250 const auto normal = crossProduct(ab, ac);
251
252 // clear list of coplanar points w.r.t to triangle ijk
253 coplanarPointBuffer.clear();
254
255 // check if this triangle ijk is admissible which means
256 // it is on the convex hull (all other points in the cloud are in the same half-space/side)
257 const bool isAdmissible = [&]()
258 {
259 // check if points are colinear and we can't form a triangle
260 // if so, skip this triangle
261 if (normal.two_norm2() < eps2*eps2)
262 return false;
263
264 int marker = 0; // 0 means undecided side (or coplanar)
265 for (int m = 0; m < numPoints; ++m)
266 {
267 if (m != i && m != j && m != k)
268 {
269 // check scalar product with surface normal to decide side
270 const auto ad = points[m] - pointI;
271 const auto sp = normal*ad;
272
273 // if the sign changes wrt the previous sign, the triangle is not part of the convex hull
274 using std::abs; using std::signbit;
275 const bool coplanar = abs(sp) < eps2*magnitude;
276 int newMarker = coplanar ? 0 : signbit(sp) ? -1 : 1;
277
278 // make decision for a side as soon as the next marker is != 0
279 // keep track of the previous marker
280 if (marker == 0 && newMarker != 0)
281 marker = newMarker;
282
283 // if marker flips, not all points are on one side
284 // zero marker (undecided side / coplanar point) shouldn't abort the process
285 if (newMarker != 0 && marker != newMarker)
286 return false;
287
288 // handle possible coplanar points
289 if (coplanar)
290 {
291 using std::abs;
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())
296 {
297 // this cluster has already been handled
298 coplanarPointBuffer.clear();
299 return false;
300 }
301 else
302 coplanarPointBuffer.push_back(points[m]);
303 }
304 }
305 }
306
307 // if there are coplanar points complete the cluster with i, j, and k
308 // and store the cluster information for future lookup
309 if (!coplanarPointBuffer.empty())
310 {
311 coplanarPointBuffer.insert(coplanarPointBuffer.end(), { points[i], points[j], points[k] });
312 coplanarClusters.emplace_back(std::make_pair(i, normal));
313 }
314
315 // we require that not all points are coplanar, so
316 // there will be always at least one non-coplanar point
317 // to check on which side the other points are.
318 // Hence, once we get here, the triangle or coplanar point cluster is part of the convex hull
319 return true;
320 }();
321
322 if (isAdmissible)
323 {
324 // check if we have a cluster of coplanar points forming on of the
325 // faces of the convex hull, if yes, compute (2d) convex hull first and triangulate
326 if (!coplanarPointBuffer.empty())
327 {
328 const auto triangles = triangulate<2, 3, TriangulationPolicy::ConvexHullPolicy>(coplanarPointBuffer);
329 for (const auto& triangle : triangles)
330 {
331 const auto ab = triangle[1] - triangle[0];
332 const auto ac = triangle[2] - triangle[0];
333 const auto normal = crossProduct(ab, ac);
334 const auto am = midPoint - triangle[0];
335 const auto sp = normal*am;
336 using std::signbit;
337 const bool isBelow = signbit(sp);
338 if (isBelow)
339 triangulation.emplace_back(Tetrahedron{
340 triangle[0], triangle[2], triangle[1], midPoint
341 });
342 else
343 triangulation.emplace_back(Tetrahedron{
344 triangle[0], triangle[1], triangle[2], midPoint
345 });
346 }
347 }
348 else
349 {
350 const auto am = midPoint - pointI;
351 const auto sp = normal*am;
352 using std::signbit;
353 const bool isBelow = signbit(sp);
354 if (isBelow)
355 triangulation.emplace_back(Tetrahedron{
356 pointI, points[k], points[j], midPoint
357 });
358 else
359 triangulation.emplace_back(Tetrahedron{
360 pointI, points[j], points[k], midPoint
361 });
362 }
363 }
364 }
365 }
366 }
367
368 // sanity check: if points are not coplanar, then using the mid point policy, we get at least 4 tetrahedrons
369 if (triangulation.size() < 4)
370 DUNE_THROW(Dune::InvalidStateException, "Something went wrong with the triangulation!");
371
372 return triangulation;
373}
374
375} // end namespace Dumux
376
377#endif
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: adapt.hh:17
Definition: triangulation.hh:39
Delaunay-type triangulations.
Definition: triangulation.hh:42
Definition: triangulation.hh:35