23#ifndef DUMUX_GRAHAM_CONVEX_HULL_HH
24#define DUMUX_GRAHAM_CONVEX_HULL_HH
31#include <dune/common/exceptions.hh>
32#include <dune/common/fvector.hh>
48 const Dune::FieldVector<ctype, 3>& b,
49 const Dune::FieldVector<ctype, 3>& c,
50 const Dune::FieldVector<ctype, 3>& normal)
55 const auto area = f*normal;
66template<
int dim,
class ctype,
67 std::enable_if_t<(dim==2),
int> = 0>
68std::vector<Dune::FieldVector<ctype, 3>>
71 using Point = Dune::FieldVector<ctype, 3>;
72 std::vector<Point> convexHull;
75 if (points.size() < 3)
79 if (points.size() == 3)
83 const auto a = points[1] - points[0];
84 auto b = points[2] - points[0];
89 auto norm = normal.two_norm();
90 while (norm == 0.0 && k < points.size()-1)
94 norm = normal.two_norm();
102 const auto eps = 1e-7*sqrt(norm);
106 auto minIt = std::min_element(points.begin(), points.end(), [&eps](
const auto& a,
const auto& b)
109 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : (abs(a[1]-b[1]) > eps ? a[1] < b[1] : (a[2] < b[2])));
113 std::iter_swap(minIt, points.begin());
117 const auto pivot = points[0];
118 std::sort(points.begin()+1, points.end(), [&](
const auto& a,
const auto& b)
120 const int order = getOrientation(pivot, a, b, normal);
122 return (a-pivot).two_norm() < (b-pivot).two_norm();
124 return (order == -1);
128 convexHull.reserve(50);
129 convexHull.push_back(points[0]);
130 convexHull.push_back(points[1]);
131 convexHull.push_back(points[2]);
136 for (std::size_t i = 3; i < points.size(); ++i)
138 Point p = convexHull.back();
139 convexHull.pop_back();
141 while (
getOrientation(convexHull.back(), p, points[i], normal) != -1)
144 if (convexHull.size() == 1)
148 assert(i < points.size()-1);
153 p = convexHull.back();
154 convexHull.pop_back();
159 convexHull.emplace_back(std::move(p));
160 convexHull.push_back(points[i]);
172template<
int dim,
class ctype,
173 std::enable_if_t<(dim==2),
int> = 0>
174std::vector<Dune::FieldVector<ctype, 2>>
177 std::vector<Dune::FieldVector<ctype, 3>> points3D;
178 points3D.reserve(points.size());
179 std::transform(points.begin(), points.end(), std::back_inserter(points3D),
180 [](
const auto& p) { return Dune::FieldVector<ctype, 3>({p[0], p[1], 0.0}); });
182 const auto result3D = grahamConvexHullImpl<2>(points3D);
184 std::vector<Dune::FieldVector<ctype, 2>> result2D;
185 result2D.reserve(result3D.size());
186 std::transform(result3D.begin(), result3D.end(), std::back_inserter(result2D),
187 [](
const auto& p) { return Dune::FieldVector<ctype, 2>({p[0], p[1]}); });
197template<
int dim,
class ctype,
int dimWorld>
198std::vector<Dune::FieldVector<ctype, dimWorld>>
grahamConvexHull(std::vector<Dune::FieldVector<ctype, dimWorld>>& points)
200 return grahamConvexHullImpl<dim>(points);
209template<
int dim,
class ctype,
int dimWorld>
210std::vector<Dune::FieldVector<ctype, dimWorld>>
grahamConvexHull(
const std::vector<Dune::FieldVector<ctype, dimWorld>>& points)
212 auto copyPoints = points;
213 return grahamConvexHullImpl<dim>(copyPoints);
Functionality to triangulate point clouds.
Define some often used mathematical functions.
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:631
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:618
std::vector< Dune::FieldVector< ctype, 3 > > grahamConvexHullImpl(std::vector< Dune::FieldVector< ctype, 3 > > &points)
Compute the points making up the convex hull around the given set of unordered points.
Definition: grahamconvexhull.hh:69
int getOrientation(const Dune::FieldVector< ctype, 3 > &a, const Dune::FieldVector< ctype, 3 > &b, const Dune::FieldVector< ctype, 3 > &c, const Dune::FieldVector< ctype, 3 > &normal)
Returns the orientation of a sequence a-->b-->c in one plane (defined by normal vector)
Definition: grahamconvexhull.hh:47
std::vector< Dune::FieldVector< ctype, dimWorld > > grahamConvexHull(const std::vector< Dune::FieldVector< ctype, dimWorld > > &points)
Compute the points making up the convex hull around the given set of unordered points.
Definition: grahamconvexhull.hh:210