3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
makegeometry.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_MAKE_GEOMETRY_HH
23#define DUMUX_MAKE_GEOMETRY_HH
24
25#include <vector>
26#include <array>
27#include <limits>
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
30#include <dune/common/exceptions.hh>
31#include <dune/geometry/multilineargeometry.hh>
32#include <dumux/common/math.hh>
34
35namespace Dumux {
36
41template<class CoordScalar>
42bool pointsAreCoplanar(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points, const CoordScalar scale)
43{
44 if (points.size() != 4)
45 DUNE_THROW(Dune::InvalidStateException, "Check only works for 4 points!");
46
47 // (see "Real-Time Collision Detection" by Christer Ericson)
48 Dune::FieldMatrix<CoordScalar, 4, 4> M;
49 for (int i = 0; i < 3; ++i )
50 M[i] = {points[0][i], points[1][i], points[2][i], points[3][i]};
51 M[3] = {1.0*scale, 1.0*scale, 1.0*scale, 1.0*scale};
52
53 using std::abs;
54 return abs(M.determinant()) < 1.5e-7*scale*scale*scale*scale;
55}
56
61template<class CoordScalar>
62bool pointsAreCoplanar(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points)
63{
64 Dune::FieldVector<CoordScalar, 3> bBoxMin(std::numeric_limits<CoordScalar>::max());
65 Dune::FieldVector<CoordScalar, 3> bBoxMax(std::numeric_limits<CoordScalar>::lowest());
66 for (const auto& p : points)
67 {
68 for (int i=0; i<3; i++)
69 {
70 using std::min;
71 using std::max;
72 bBoxMin[i] = min(bBoxMin[i], p[i]);
73 bBoxMax[i] = max(bBoxMax[i], p[i]);
74 }
75 }
76
77 const auto size = (bBoxMax - bBoxMin).two_norm();
78
79 return pointsAreCoplanar(points, size);
80}
81
89template<class CoordScalar>
90std::vector<Dune::FieldVector<CoordScalar, 3>> getReorderedPoints(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points)
91{
92 std::array<int, 4> tmp;
93 return getReorderedPoints(points, tmp);
94}
95
103template<class CoordScalar>
104std::vector<Dune::FieldVector<CoordScalar, 3>> getReorderedPoints(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points,
105 std::array<int, 4>& orientations)
106{
107 if(points.size() == 4)
108 {
109 auto& p0 = points[0];
110 auto& p1 = points[1];
111 auto& p2 = points[2];
112 auto& p3 = points[3];
113
114 // check if the points define a proper quadrilateral
115 const auto normal = crossProduct((p1 - p0), (p2 - p0));
116
117 orientations = { getOrientation(p0, p3, p2, normal),
118 getOrientation(p0, p3, p1, normal),
119 getOrientation(p2, p1, p0, normal),
120 getOrientation(p2, p1, p3, normal) };
121
122
123 // check if the points follow the dune ordering (see http://www.dcs.gla.ac.uk/~pat/52233/slides/Geometry1x1.pdf)
124 const bool diagonalsIntersect = (orientations[0] != orientations[1]) && (orientations[2] != orientations[3]);
125
126 // the points conform with the dune ordering
127 if(diagonalsIntersect)
128 return points;
129
130 // the points do not conform with the dune ordering, re-order
131 using GlobalPosition = Dune::FieldVector<CoordScalar, 3>;
132 if(!diagonalsIntersect && orientations[0] == 1)
133 return std::vector<GlobalPosition>{p1, p0, p2, p3};
134 else if(!diagonalsIntersect && orientations[0] == -1)
135 return std::vector<GlobalPosition>{p3, p1, p0, p2};
136 else
137 DUNE_THROW(Dune::InvalidStateException, "Could not reorder points");
138 }
139 else
140 DUNE_THROW(Dune::NotImplemented, "Reorder for " << points.size() << " points.");
141}
142
151template<class CoordScalar, bool enableSanityCheck = true>
152auto makeDuneQuadrilaterial(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points)
153{
154 if (points.size() != 4)
155 DUNE_THROW(Dune::InvalidStateException, "A quadrilateral needs 4 corner points!");
156
157 using GlobalPosition = Dune::FieldVector<CoordScalar, 3>;
158 static constexpr auto coordDim = GlobalPosition::dimension;
159 static constexpr auto dim = coordDim-1;
160 using GeometryType = Dune::MultiLinearGeometry<CoordScalar, dim, coordDim>;
161
162 // if no sanity check if desired, use the given points directly to construct the geometry
163 if (!enableSanityCheck)
164 return GeometryType(Dune::GeometryTypes::quadrilateral, points);
165
166 // compute size
167 Dune::FieldVector<CoordScalar, 3> bBoxMin(std::numeric_limits<CoordScalar>::max());
168 Dune::FieldVector<CoordScalar, 3> bBoxMax(std::numeric_limits<CoordScalar>::lowest());
169 for (const auto& p : points)
170 {
171 for (int i = 0; i < 3; i++)
172 {
173 using std::min;
174 using std::max;
175 bBoxMin[i] = min(bBoxMin[i], p[i]);
176 bBoxMax[i] = max(bBoxMax[i], p[i]);
177 }
178 }
179
180 const auto size = (bBoxMax - bBoxMin).two_norm();
181
182 // otherwise, perform a number of checks and corrections
183 if (!pointsAreCoplanar(points, size))
184 DUNE_THROW(Dune::InvalidStateException, "Points do not lie within a plane");
185
186 auto corners = grahamConvexHull<2>(points);
187 if (corners.size() != 4)
188 DUNE_THROW(Dune::InvalidStateException, "Points do not span a strictly convex polygon!");
189
190 // make sure points conform with dune ordering
191 std::array<int, 4> orientations;
192 corners = getReorderedPoints(corners, orientations);
193
194 if (std::any_of(orientations.begin(), orientations.end(), [](auto i){ return i == 0; }))
195 DUNE_THROW(Dune::InvalidStateException, "More than two points lie on the same line.");
196
197 const auto quadrilateral = GeometryType(Dune::GeometryTypes::quadrilateral, corners);
198
199 const auto eps = 1e-7;
200 if (quadrilateral.volume() < eps*size*size)
201 DUNE_THROW(Dune::InvalidStateException, "Something went wrong, geometry has area of zero");
202
203 return quadrilateral;
204}
205
206
207
208} // end namespace Dumux
209
210#endif // DUMUX_MAKE_GEOMETRY_HH
A function to compute the convex hull of a point cloud and a function to triangulate the polygon area...
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
auto makeDuneQuadrilaterial(const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points)
Creates a dune quadrilateral geometry given 4 corner points.
Definition: makegeometry.hh:152
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< CoordScalar, 3 > > getReorderedPoints(const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points)
Returns a vector of points following the dune ordering. Convenience method that creates a temporary o...
Definition: makegeometry.hh:90
bool pointsAreCoplanar(const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points, const CoordScalar scale)
Checks if four points lie within the same plane.
Definition: makegeometry.hh:42
Definition: adapt.hh:29