3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
grahamconvexhull.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 *****************************************************************************/
23#ifndef DUMUX_GEOMETRY_GRAHAM_CONVEX_HULL_HH
24#define DUMUX_GEOMETRY_GRAHAM_CONVEX_HULL_HH
25
26#include <vector>
27#include <array>
28#include <algorithm>
29#include <iterator>
30
31#include <dune/common/exceptions.hh>
32#include <dune/common/fvector.hh>
33
34#include <dumux/common/math.hh>
36
37namespace Dumux {
38
46template<class ctype>
47int getOrientation(const Dune::FieldVector<ctype, 3>& a,
48 const Dune::FieldVector<ctype, 3>& b,
49 const Dune::FieldVector<ctype, 3>& c,
50 const Dune::FieldVector<ctype, 3>& normal)
51{
52 const auto d = b-a;
53 const auto e = c-b;
54 const auto f = Dumux::crossProduct(d, e);
55 const auto area = f*normal;
56 return Dumux::sign(-area);
57}
58
66template<int dim, class ctype,
67 std::enable_if_t<(dim==2), int> = 0>
68std::vector<Dune::FieldVector<ctype, 3>>
69grahamConvexHullImpl(std::vector<Dune::FieldVector<ctype, 3>>& points)
70{
71 using Point = Dune::FieldVector<ctype, 3>;
72 std::vector<Point> convexHull;
73
74 // return empty convex hull
75 if (points.size() < 3)
76 return convexHull;
77
78 // return the points (already just one triangle)
79 if (points.size() == 3)
80 return points;
81
82 // try to compute the normal vector of the plane
83 const auto a = points[1] - points[0];
84 auto b = points[2] - points[0];
85 auto normal = Dumux::crossProduct(a, b);
86
87 // make sure the normal vector has non-zero length
88 std::size_t k = 2;
89 auto norm = normal.two_norm();
90 while (norm == 0.0 && k < points.size()-1)
91 {
92 b = points[++k];
94 norm = normal.two_norm();
95 }
96
97 // if all given points are colinear -> return empty convex hull
98 if (norm == 0.0)
99 return convexHull;
100
101 using std::sqrt;
102 const auto eps = 1e-7*sqrt(norm);
103 normal /= norm;
104
105 // find the element with the smallest x coordinate (if x is the same, smallest y coordinate, and so on...)
106 auto minIt = std::min_element(points.begin(), points.end(), [&eps](const auto& a, const auto& b)
107 {
108 using std::abs;
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])));
110 });
111
112 // swap the smallest element to the front
113 std::iter_swap(minIt, points.begin());
114
115 // choose the first (min element) as the pivot point
116 // sort in counter-clockwise order around pivot point
117 const auto pivot = points[0];
118 std::sort(points.begin()+1, points.end(), [&](const auto& a, const auto& b)
119 {
120 const int order = getOrientation(pivot, a, b, normal);
121 if (order == 0)
122 return (a-pivot).two_norm() < (b-pivot).two_norm();
123 else
124 return (order == -1);
125 });
126
127 // push the first three points
128 convexHull.reserve(50);
129 convexHull.push_back(points[0]);
130 convexHull.push_back(points[1]);
131 convexHull.push_back(points[2]);
132
133 // This is the heart of the algorithm
134 // pop_back until the last point in the queue forms a counter-clockwise oriented line
135 // with the next two vertices. Then add these points to the queue.
136 for (std::size_t i = 3; i < points.size(); ++i)
137 {
138 Point p = convexHull.back();
139 convexHull.pop_back();
140 // keep popping until the orientation a->b->currentp is counter-clockwise
141 while (getOrientation(convexHull.back(), p, points[i], normal) != -1)
142 {
143 // make sure the queue doesn't get empty
144 if (convexHull.size() == 1)
145 {
146 // before we reach i=size-1 there has to be a good candidate
147 // as not all points are colinear (a non-zero plane normal exists)
148 assert(i < points.size()-1);
149 p = points[i++];
150 }
151 else
152 {
153 p = convexHull.back();
154 convexHull.pop_back();
155 }
156 }
157
158 // add back the last popped point and this point
159 convexHull.emplace_back(std::move(p));
160 convexHull.push_back(points[i]);
161 }
162
163 return convexHull;
164}
165
172template<int dim, class ctype,
173 std::enable_if_t<(dim==2), int> = 0>
174std::vector<Dune::FieldVector<ctype, 2>>
175grahamConvexHullImpl(const std::vector<Dune::FieldVector<ctype, 2>>& points)
176{
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}); });
181
182 const auto result3D = grahamConvexHullImpl<2>(points3D);
183
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]}); });
188
189 return result2D;
190}
191
197template<int dim, class ctype, int dimWorld>
198std::vector<Dune::FieldVector<ctype, dimWorld>> grahamConvexHull(std::vector<Dune::FieldVector<ctype, dimWorld>>& points)
199{
200 return grahamConvexHullImpl<dim>(points);
201}
202
209template<int dim, class ctype, int dimWorld>
210std::vector<Dune::FieldVector<ctype, dimWorld>> grahamConvexHull(const std::vector<Dune::FieldVector<ctype, dimWorld>>& points)
211{
212 auto copyPoints = points;
213 return grahamConvexHullImpl<dim>(copyPoints);
214}
215
216} // end namespace Dumux
217
218#endif
Define some often used mathematical functions.
Functionality to triangulate point clouds.
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
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< 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
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
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition: math.hh:641
Definition: adapt.hh:29