3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
makegeometry.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3/*****************************************************************************
4 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_GEOMETRY_MAKE_GEOMETRY_HH
25#define DUMUX_GEOMETRY_MAKE_GEOMETRY_HH
26
27#include <vector>
28#include <array>
29#include <limits>
30#include <dune/common/fvector.hh>
31#include <dune/common/fmatrix.hh>
32#include <dune/common/exceptions.hh>
33#include <dune/geometry/multilineargeometry.hh>
34#include <dumux/common/math.hh>
36
37namespace Dumux {
38
43template<class CoordScalar>
44bool pointsAreCoplanar(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points, const CoordScalar scale)
45{
46 if (points.size() != 4)
47 DUNE_THROW(Dune::InvalidStateException, "Check only works for 4 points!");
48
49 // (see "Real-Time Collision Detection" by Christer Ericson)
50 Dune::FieldMatrix<CoordScalar, 4, 4> M;
51 for (int i = 0; i < 3; ++i )
52 M[i] = {points[0][i], points[1][i], points[2][i], points[3][i]};
53 M[3] = {1.0*scale, 1.0*scale, 1.0*scale, 1.0*scale};
54
55 using std::abs;
56 return abs(M.determinant()) < 1.5e-7*scale*scale*scale*scale;
57}
58
63template<class CoordScalar>
64bool pointsAreCoplanar(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points)
65{
66 Dune::FieldVector<CoordScalar, 3> bBoxMin(std::numeric_limits<CoordScalar>::max());
67 Dune::FieldVector<CoordScalar, 3> bBoxMax(std::numeric_limits<CoordScalar>::lowest());
68 for (const auto& p : points)
69 {
70 for (int i=0; i<3; i++)
71 {
72 using std::min;
73 using std::max;
74 bBoxMin[i] = min(bBoxMin[i], p[i]);
75 bBoxMax[i] = max(bBoxMax[i], p[i]);
76 }
77 }
78
79 const auto size = (bBoxMax - bBoxMin).two_norm();
80
81 return pointsAreCoplanar(points, size);
82}
83
91template<class CoordScalar>
92std::vector<Dune::FieldVector<CoordScalar, 3>> getReorderedPoints(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points)
93{
94 std::array<int, 4> tmp;
95 return getReorderedPoints(points, tmp);
96}
97
105template<class CoordScalar>
106std::vector<Dune::FieldVector<CoordScalar, 3>> getReorderedPoints(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points,
107 std::array<int, 4>& orientations)
108{
109 if(points.size() == 4)
110 {
111 auto& p0 = points[0];
112 auto& p1 = points[1];
113 auto& p2 = points[2];
114 auto& p3 = points[3];
115
116 // check if the points define a proper quadrilateral
117 const auto normal = crossProduct((p1 - p0), (p2 - p0));
118
119 orientations = { getOrientation(p0, p3, p2, normal),
120 getOrientation(p0, p3, p1, normal),
121 getOrientation(p2, p1, p0, normal),
122 getOrientation(p2, p1, p3, normal) };
123
124
125 // check if the points follow the dune ordering (see http://www.dcs.gla.ac.uk/~pat/52233/slides/Geometry1x1.pdf)
126 const bool diagonalsIntersect = (orientations[0] != orientations[1]) && (orientations[2] != orientations[3]);
127
128 // the points conform with the dune ordering
129 if(diagonalsIntersect)
130 return points;
131
132 // the points do not conform with the dune ordering, re-order
133 using GlobalPosition = Dune::FieldVector<CoordScalar, 3>;
134 if(!diagonalsIntersect && orientations[0] == 1)
135 return std::vector<GlobalPosition>{p1, p0, p2, p3};
136 else if(!diagonalsIntersect && orientations[0] == -1)
137 return std::vector<GlobalPosition>{p3, p1, p0, p2};
138 else
139 DUNE_THROW(Dune::InvalidStateException, "Could not reorder points");
140 }
141 else
142 DUNE_THROW(Dune::NotImplemented, "Reorder for " << points.size() << " points.");
143}
144
153template<class CoordScalar, bool enableSanityCheck = true>
154auto makeDuneQuadrilaterial(const std::vector<Dune::FieldVector<CoordScalar, 3>>& points)
155{
156 if (points.size() != 4)
157 DUNE_THROW(Dune::InvalidStateException, "A quadrilateral needs 4 corner points!");
158
159 using GlobalPosition = Dune::FieldVector<CoordScalar, 3>;
160 static constexpr auto coordDim = GlobalPosition::dimension;
161 static constexpr auto dim = coordDim-1;
162 using GeometryType = Dune::MultiLinearGeometry<CoordScalar, dim, coordDim>;
163
164 // if no sanity check if desired, use the given points directly to construct the geometry
165 if (!enableSanityCheck)
166 return GeometryType(Dune::GeometryTypes::quadrilateral, points);
167
168 // compute size
169 Dune::FieldVector<CoordScalar, 3> bBoxMin(std::numeric_limits<CoordScalar>::max());
170 Dune::FieldVector<CoordScalar, 3> bBoxMax(std::numeric_limits<CoordScalar>::lowest());
171 for (const auto& p : points)
172 {
173 for (int i = 0; i < 3; i++)
174 {
175 using std::min;
176 using std::max;
177 bBoxMin[i] = min(bBoxMin[i], p[i]);
178 bBoxMax[i] = max(bBoxMax[i], p[i]);
179 }
180 }
181
182 const auto size = (bBoxMax - bBoxMin).two_norm();
183
184 // otherwise, perform a number of checks and corrections
185 if (!pointsAreCoplanar(points, size))
186 DUNE_THROW(Dune::InvalidStateException, "Points do not lie within a plane");
187
188 auto corners = grahamConvexHull<2>(points);
189 if (corners.size() != 4)
190 DUNE_THROW(Dune::InvalidStateException, "Points do not span a strictly convex polygon!");
191
192 // make sure points conform with dune ordering
193 std::array<int, 4> orientations;
194 corners = getReorderedPoints(corners, orientations);
195
196 if (std::any_of(orientations.begin(), orientations.end(), [](auto i){ return i == 0; }))
197 DUNE_THROW(Dune::InvalidStateException, "More than two points lie on the same line.");
198
199 const auto quadrilateral = GeometryType(Dune::GeometryTypes::quadrilateral, corners);
200
201 const auto eps = 1e-7;
202 if (quadrilateral.volume() < eps*size*size)
203 DUNE_THROW(Dune::InvalidStateException, "Something went wrong, geometry has area of zero");
204
205 return quadrilateral;
206}
207
208
209
210} // end namespace Dumux
211
212#endif
Define some often used mathematical functions.
A function to compute the convex hull of a point cloud.
auto makeDuneQuadrilaterial(const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points)
Creates a dune quadrilateral geometry given 4 corner points.
Definition: makegeometry.hh:154
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:38
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:92
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:44
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29