3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
triangulation.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
24#ifndef DUMUX_GEOMETRY_TRIANGULATION_HH
25#define DUMUX_GEOMETRY_TRIANGULATION_HH
26
27#include <vector>
28#include <array>
29#include <algorithm>
30#include <type_traits>
31#include <tuple>
32
33#include <dune/common/exceptions.hh>
34#include <dune/common/fvector.hh>
35
36#include <dumux/common/math.hh>
39
40namespace Dumux {
41namespace TriangulationPolicy {
42
46
50
53
54#ifndef DOXYGEN
55namespace Detail {
56using DefaultDimPolicies = std::tuple<DelaunayPolicy, MidPointPolicy, ConvexHullPolicy>;
57} // end namespace Detail
58#endif
59
61template<int dim, int dimWorld>
62using DefaultPolicy = std::tuple_element_t<dim-1, Detail::DefaultDimPolicies>;
63
64} // end namespace TriangulationPolicy
65
73template<int dim, int dimWorld, class ctype>
74using Triangulation = std::vector< std::array<Dune::FieldVector<ctype, dimWorld>, dim+1> >;
75
86template< int dim, int dimWorld, class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>,
87 class RandomAccessContainer,
88 std::enable_if_t< std::is_same_v<Policy, TriangulationPolicy::DelaunayPolicy>
89 && dim == 1, int> = 0 >
91triangulate(const RandomAccessContainer& points)
92{
93 using ctype = typename RandomAccessContainer::value_type::value_type;
94 using Point = Dune::FieldVector<ctype, dimWorld>;
95
96 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
97 "Triangulation expects Dune::FieldVector as point type");
98
99 if (points.size() == 2)
100 return Triangulation<dim, dimWorld, ctype>({ {points[0], points[1]} });
101
103 assert(points.size() > 1);
104 DUNE_THROW(Dune::NotImplemented, "1d triangulation for point cloud size > 2");
105}
106
119template< int dim, int dimWorld, class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>,
120 class RandomAccessContainer,
121 std::enable_if_t< std::is_same_v<Policy, TriangulationPolicy::MidPointPolicy>
122 && dim == 2, int> = 0 >
123inline Triangulation<dim, dimWorld, typename RandomAccessContainer::value_type::value_type>
124triangulate(const RandomAccessContainer& convexHullPoints)
125{
126 using ctype = typename RandomAccessContainer::value_type::value_type;
127 using Point = Dune::FieldVector<ctype, dimWorld>;
128 using Triangle = std::array<Point, 3>;
129
130 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
131 "Triangulation expects Dune::FieldVector as point type");
132
133 if (convexHullPoints.size() < 3)
134 DUNE_THROW(Dune::InvalidStateException, "Try to triangulate point cloud with less than 3 points!");
135
136 if (convexHullPoints.size() == 3)
137 return std::vector<Triangle>(1, {convexHullPoints[0], convexHullPoints[1], convexHullPoints[2]});
138
139 Point midPoint(0.0);
140 for (const auto& p : convexHullPoints)
141 midPoint += p;
142 midPoint /= convexHullPoints.size();
143
144 std::vector<Triangle> triangulation;
145 triangulation.reserve(convexHullPoints.size());
146
147 for (std::size_t i = 0; i < convexHullPoints.size()-1; ++i)
148 triangulation.emplace_back(Triangle{midPoint, convexHullPoints[i], convexHullPoints[i+1]});
149
150 triangulation.emplace_back(Triangle{midPoint, convexHullPoints[convexHullPoints.size()-1], convexHullPoints[0]});
151
152 return triangulation;
153}
154
167template< int dim, int dimWorld, class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>,
168 class RandomAccessContainer,
169 std::enable_if_t< std::is_same_v<Policy, TriangulationPolicy::ConvexHullPolicy>
170 && dim == 2, int> = 0 >
171inline Triangulation<dim, dimWorld, typename RandomAccessContainer::value_type::value_type>
172triangulate(const RandomAccessContainer& points)
173{
174 const auto convexHullPoints = grahamConvexHull<2>(points);
175 return triangulate<dim, dimWorld, TriangulationPolicy::MidPointPolicy>(convexHullPoints);
176}
177
188template< int dim, int dimWorld, class Policy = TriangulationPolicy::DefaultPolicy<dim, dimWorld>,
189 class RandomAccessContainer,
190 std::enable_if_t< std::is_same_v<Policy, TriangulationPolicy::ConvexHullPolicy>
191 && dim == 3, int> = 0 >
192inline Triangulation<dim, dimWorld, typename RandomAccessContainer::value_type::value_type>
193triangulate(const RandomAccessContainer& points)
194{
195 using ctype = typename RandomAccessContainer::value_type::value_type;
196 using Point = Dune::FieldVector<ctype, dimWorld>;
197 using Tetrahedron = std::array<Point, 4>;
198
199 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
200 "Triangulation expects Dune::FieldVector as point type");
201
202 const auto numPoints = points.size();
203 if (numPoints < 4)
204 DUNE_THROW(Dune::InvalidStateException, "Trying to create 3D triangulation of point cloud with less than 4 points!");
205
206 if (numPoints == 4)
207 return std::vector<Tetrahedron>(1, {points[0], points[1], points[2], points[3]});
208
209 // compute the mid point of the point cloud (not the midpoint of the convex hull but this
210 // should not matter too much for the applications we have in mind here)
211 Point midPoint(0.0);
212 Point lowerLeft(1e100);
213 Point upperRight(-1e100);
214 for (const auto& p : points)
215 {
216 midPoint += p;
217 for (int i = 0; i < dimWorld; ++i)
218 {
219 using std::max; using std::min;
220 lowerLeft[i] = min(p[i], lowerLeft[i]);
221 upperRight[i] = max(p[i], upperRight[i]);
222 }
223 }
224 midPoint /= numPoints;
225
226 auto magnitude = 0.0;
227 using std::max;
228 for (int i = 0; i < dimWorld; ++i)
229 magnitude = max(upperRight[i] - lowerLeft[i], magnitude);
230 const auto eps = 1e-7*magnitude;
231 const auto eps2 = eps*magnitude;
232
233 // reserve memory conservatively to avoid reallocation
234 std::vector<Tetrahedron> triangulation;
235 triangulation.reserve(numPoints);
236
237 // make a buffer for storing coplanar points and indices
238 std::vector<Point> coplanarPointBuffer;
239 coplanarPointBuffer.reserve(std::min<std::size_t>(12, numPoints-1));
240
241 // remember coplanar cluster planes
242 // coplanar clusters are uniquely identified by a point and a normal (plane)
243 // we only want to add each cluster once (when handling the first triangle in the cluster)
244 std::vector<std::pair<int, Point>> coplanarClusters;
245 coplanarClusters.reserve(numPoints/3);
246
247 // brute force algorithm: Try all possible triangles and check
248 // if they are triangles of the convex hull. This is achieved by checking if all
249 // other points are in the same half-space (side)
250 // the algorithm is O(n^4), so only do this for very small point clouds
251 for (int i = 0; i < numPoints; ++i)
252 {
253 for (int j = i+1; j < numPoints; ++j)
254 {
255 for (int k = j+1; k < numPoints; ++k)
256 {
257 const auto pointI = points[i];
258 const auto ab = points[j] - pointI;
259 const auto ac = points[k] - pointI;
260 const auto normal = crossProduct(ab, ac);
261
262 // clear list of coplanar points w.r.t to triangle ijk
263 coplanarPointBuffer.clear();
264
265 // check if this triangle ijk is admissible which means
266 // it is on the convex hull (all other points in the cloud are in the same half-space/side)
267 const bool isAdmissible = [&]()
268 {
269 // check if points are colinear and we can't form a triangle
270 // if so, skip this triangle
271 if (normal.two_norm2() < eps2*eps2)
272 return false;
273
274 int marker = 0; // 0 means undecided side (or coplanar)
275 for (int m = 0; m < numPoints; ++m)
276 {
277 if (m != i && m != j && m != k)
278 {
279 // check scalar product with surface normal to decide side
280 const auto ad = points[m] - pointI;
281 const auto sp = normal*ad;
282
283 // if the sign changes wrt the previous sign, the triangle is not part of the convex hull
284 using std::abs; using std::signbit;
285 const bool coplanar = abs(sp) < eps2*magnitude;
286 int newMarker = coplanar ? 0 : signbit(sp) ? -1 : 1;
287
288 // make decision for a side as soon as the next marker is != 0
289 // keep track of the previous marker
290 if (marker == 0 && newMarker != 0)
291 marker = newMarker;
292
293 // if marker flips, not all points are on one side
294 // zero marker (undecided side / coplanar point) shouldn't abort the process
295 if (newMarker != 0 && marker != newMarker)
296 return false;
297
298 // handle possible coplanar points
299 if (coplanar)
300 {
301 using std::abs;
302 if (m < k && std::find_if(
303 coplanarClusters.begin(), coplanarClusters.end(),
304 [=](const auto& c){ return c.first == std::min(m, i) && abs(ab*c.second) < eps2*magnitude; }
305 ) != coplanarClusters.end())
306 {
307 // this cluster has already been handled
308 coplanarPointBuffer.clear();
309 return false;
310 }
311 else
312 coplanarPointBuffer.push_back(points[m]);
313 }
314 }
315 }
316
317 // if there are coplanar points complete the cluster with i, j, and k
318 // and store the cluster information for future lookup
319 if (!coplanarPointBuffer.empty())
320 {
321 coplanarPointBuffer.insert(coplanarPointBuffer.end(), { points[i], points[j], points[k] });
322 coplanarClusters.emplace_back(std::make_pair(i, normal));
323 }
324
325 // we require that not all points are coplanar, so
326 // there will be always at least one non-coplanar point
327 // to check on which side the other points are.
328 // Hence, once we get here, the triangle or coplanar point cluster is part of the convex hull
329 return true;
330 }();
331
332 if (isAdmissible)
333 {
334 // check if we have a cluster of coplanar points forming on of the
335 // faces of the convex hull, if yes, compute (2d) convex hull first and triangulate
336 if (!coplanarPointBuffer.empty())
337 {
338 const auto triangles = triangulate<2, 3, TriangulationPolicy::ConvexHullPolicy>(coplanarPointBuffer);
339 for (const auto& triangle : triangles)
340 {
341 const auto ab = triangle[1] - triangle[0];
342 const auto ac = triangle[2] - triangle[0];
343 const auto normal = crossProduct(ab, ac);
344 const auto am = midPoint - triangle[0];
345 const auto sp = normal*am;
346 using std::signbit;
347 const bool isBelow = signbit(sp);
348 if (isBelow)
349 triangulation.emplace_back(Tetrahedron{
350 triangle[0], triangle[2], triangle[1], midPoint
351 });
352 else
353 triangulation.emplace_back(Tetrahedron{
354 triangle[0], triangle[1], triangle[2], midPoint
355 });
356 }
357 }
358 else
359 {
360 const auto am = midPoint - pointI;
361 const auto sp = normal*am;
362 using std::signbit;
363 const bool isBelow = signbit(sp);
364 if (isBelow)
365 triangulation.emplace_back(Tetrahedron{
366 pointI, points[k], points[j], midPoint
367 });
368 else
369 triangulation.emplace_back(Tetrahedron{
370 pointI, points[j], points[k], midPoint
371 });
372 }
373 }
374 }
375 }
376 }
377
378 // sanity check: if points are not coplanar, then using the mid point policy, we get at least 4 tetrahedrons
379 if (triangulation.size() < 4)
380 DUNE_THROW(Dune::InvalidStateException, "Something went wrong with the triangulation!");
381
382 return triangulation;
383}
384
385} // end namespace Dumux
386
387#endif
Define some often used mathematical functions.
Detect if a point intersects a simplex (including boundary)
A function to compute the convex hull of a point cloud.
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:91
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:36
std::vector< std::array< Dune::FieldVector< ctype, dimWorld >, dim+1 > > Triangulation
The default data type to store triangulations.
Definition: triangulation.hh:74
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:654
Definition: adapt.hh:29
std::tuple_element_t< dim-1, Detail::DefaultDimPolicies > DefaultPolicy
Default policy for a given dimension.
Definition: triangulation.hh:62
Definition: triangulation.hh:45
Definition: triangulation.hh:49
Delaunay-type triangulations.
Definition: triangulation.hh:52