24#ifndef DUMUX_GEOMETRY_TRIANGULATION_HH
25#define DUMUX_GEOMETRY_TRIANGULATION_HH
33#include <dune/common/exceptions.hh>
34#include <dune/common/fvector.hh>
41namespace TriangulationPolicy {
56using DefaultDimPolicies = std::tuple<DelaunayPolicy, MidPointPolicy, ConvexHullPolicy>;
61template<
int dim,
int dimWorld>
62using DefaultPolicy = std::tuple_element_t<dim-1, Detail::DefaultDimPolicies>;
73template<
int dim,
int dimWorld,
class ctype>
74using Triangulation = std::vector< std::array<Dune::FieldVector<ctype, dimWorld>, dim+1> >;
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 >
93 using ctype =
typename RandomAccessContainer::value_type::value_type;
94 using Point = Dune::FieldVector<ctype, dimWorld>;
96 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
97 "Triangulation expects Dune::FieldVector as point type");
99 if (points.size() == 2)
103 assert(points.size() > 1);
104 DUNE_THROW(Dune::NotImplemented,
"1d triangulation for point cloud size > 2");
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::M
idPo
intPolicy>
122 && dim == 2,
int> = 0 >
123inline Triangulation<dim, dimWorld, typename RandomAccessContainer::value_type::value_type>
124triangulate(
const RandomAccessContainer& convexHullPoints)
126 using ctype =
typename RandomAccessContainer::value_type::value_type;
127 using Point = Dune::FieldVector<ctype, dimWorld>;
128 using Triangle = std::array<Point, 3>;
130 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
131 "Triangulation expects Dune::FieldVector as point type");
133 if (convexHullPoints.size() < 3)
134 DUNE_THROW(Dune::InvalidStateException,
"Try to triangulate point cloud with less than 3 points!");
136 if (convexHullPoints.size() == 3)
137 return std::vector<Triangle>(1, {convexHullPoints[0], convexHullPoints[1], convexHullPoints[2]});
140 for (
const auto& p : convexHullPoints)
142 midPoint /= convexHullPoints.size();
144 std::vector<Triangle> triangulation;
145 triangulation.reserve(convexHullPoints.size());
147 for (std::size_t i = 0; i < convexHullPoints.size()-1; ++i)
148 triangulation.emplace_back(Triangle{midPoint, convexHullPoints[i], convexHullPoints[i+1]});
150 triangulation.emplace_back(Triangle{midPoint, convexHullPoints[convexHullPoints.size()-1], convexHullPoints[0]});
152 return triangulation;
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>
174 const auto convexHullPoints = grahamConvexHull<2>(points);
175 return triangulate<dim, dimWorld, TriangulationPolicy::MidPointPolicy>(convexHullPoints);
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>
195 using ctype =
typename RandomAccessContainer::value_type::value_type;
196 using Point = Dune::FieldVector<ctype, dimWorld>;
197 using Tetrahedron = std::array<Point, 4>;
199 static_assert(std::is_same_v<typename RandomAccessContainer::value_type, Point>,
200 "Triangulation expects Dune::FieldVector as point type");
202 const auto numPoints = points.size();
204 DUNE_THROW(Dune::InvalidStateException,
"Trying to create 3D triangulation of point cloud with less than 4 points!");
207 return std::vector<Tetrahedron>(1, {points[0], points[1], points[2], points[3]});
212 Point lowerLeft(1e100);
213 Point upperRight(-1e100);
214 for (
const auto& p : points)
217 for (
int i = 0; i < dimWorld; ++i)
219 using std::max;
using std::min;
220 lowerLeft[i] = min(p[i], lowerLeft[i]);
221 upperRight[i] = max(p[i], upperRight[i]);
224 midPoint /= numPoints;
226 auto magnitude = 0.0;
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;
234 std::vector<Tetrahedron> triangulation;
235 triangulation.reserve(numPoints);
238 std::vector<Point> coplanarPointBuffer;
239 coplanarPointBuffer.reserve(std::min<std::size_t>(12, numPoints-1));
244 std::vector<std::pair<int, Point>> coplanarClusters;
245 coplanarClusters.reserve(numPoints/3);
251 for (
int i = 0; i < numPoints; ++i)
253 for (
int j = i+1; j < numPoints; ++j)
255 for (
int k = j+1; k < numPoints; ++k)
257 const auto pointI = points[i];
258 const auto ab = points[j] - pointI;
259 const auto ac = points[k] - pointI;
263 coplanarPointBuffer.clear();
267 const bool isAdmissible = [&]()
271 if (
normal.two_norm2() < eps2*eps2)
275 for (
int m = 0; m < numPoints; ++m)
277 if (m != i && m != j && m != k)
280 const auto ad = points[m] - pointI;
281 const auto sp =
normal*ad;
284 using std::abs;
using std::signbit;
285 const bool coplanar = abs(sp) < eps2*magnitude;
286 int newMarker = coplanar ? 0 : signbit(sp) ? -1 : 1;
290 if (marker == 0 && newMarker != 0)
295 if (newMarker != 0 && marker != newMarker)
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())
308 coplanarPointBuffer.clear();
312 coplanarPointBuffer.push_back(points[m]);
319 if (!coplanarPointBuffer.empty())
321 coplanarPointBuffer.insert(coplanarPointBuffer.end(), { points[i], points[j], points[k] });
322 coplanarClusters.emplace_back(std::make_pair(i,
normal));
336 if (!coplanarPointBuffer.empty())
338 const auto triangles = triangulate<2, 3, TriangulationPolicy::ConvexHullPolicy>(coplanarPointBuffer);
339 for (
const auto& triangle : triangles)
341 const auto ab = triangle[1] - triangle[0];
342 const auto ac = triangle[2] - triangle[0];
344 const auto am = midPoint - triangle[0];
345 const auto sp =
normal*am;
347 const bool isBelow = signbit(sp);
349 triangulation.emplace_back(Tetrahedron{
350 triangle[0], triangle[2], triangle[1], midPoint
353 triangulation.emplace_back(Tetrahedron{
354 triangle[0], triangle[1], triangle[2], midPoint
360 const auto am = midPoint - pointI;
361 const auto sp =
normal*am;
363 const bool isBelow = signbit(sp);
365 triangulation.emplace_back(Tetrahedron{
366 pointI, points[k], points[j], midPoint
369 triangulation.emplace_back(Tetrahedron{
370 pointI, points[j], points[k], midPoint
379 if (triangulation.size() < 4)
380 DUNE_THROW(Dune::InvalidStateException,
"Something went wrong with the triangulation!");
382 return triangulation;
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
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