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