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