24#ifndef DUMUX_GEOMETRY_MAKE_GEOMETRY_HH
25#define DUMUX_GEOMETRY_MAKE_GEOMETRY_HH
30#include <dune/common/fvector.hh>
31#include <dune/common/fmatrix.hh>
32#include <dune/common/exceptions.hh>
33#include <dune/geometry/multilineargeometry.hh>
43template<
class CoordScalar>
44bool pointsAreCoplanar(
const std::vector<Dune::FieldVector<CoordScalar, 3>>& points,
const CoordScalar scale)
46 if (points.size() != 4)
47 DUNE_THROW(Dune::InvalidStateException,
"Check only works for 4 points!");
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};
56 return abs(M.determinant()) < 1.5e-7*scale*scale*scale*scale;
63template<
class CoordScalar>
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)
70 for (
int i=0; i<3; i++)
74 bBoxMin[i] = min(bBoxMin[i], p[i]);
75 bBoxMax[i] = max(bBoxMax[i], p[i]);
79 const auto size = (bBoxMax - bBoxMin).two_norm();
91template<
class CoordScalar>
92std::vector<Dune::FieldVector<CoordScalar, 3>>
getReorderedPoints(
const std::vector<Dune::FieldVector<CoordScalar, 3>>& points)
94 std::array<int, 4> tmp;
105template<
class CoordScalar>
106std::vector<Dune::FieldVector<CoordScalar, 3>>
getReorderedPoints(
const std::vector<Dune::FieldVector<CoordScalar, 3>>& points,
107 std::array<int, 4>& orientations)
109 if(points.size() == 4)
111 auto& p0 = points[0];
112 auto& p1 = points[1];
113 auto& p2 = points[2];
114 auto& p3 = points[3];
126 const bool diagonalsIntersect = (orientations[0] != orientations[1]) && (orientations[2] != orientations[3]);
129 if(diagonalsIntersect)
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};
139 DUNE_THROW(Dune::InvalidStateException,
"Could not reorder points");
142 DUNE_THROW(Dune::NotImplemented,
"Reorder for " << points.size() <<
" points.");
153template<
class CoordScalar,
bool enableSanityCheck = true>
156 if (points.size() != 4)
157 DUNE_THROW(Dune::InvalidStateException,
"A quadrilateral needs 4 corner points!");
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>;
165 if (!enableSanityCheck)
166 return GeometryType(Dune::GeometryTypes::quadrilateral, points);
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)
173 for (
int i = 0; i < 3; i++)
177 bBoxMin[i] = min(bBoxMin[i], p[i]);
178 bBoxMax[i] = max(bBoxMax[i], p[i]);
182 const auto size = (bBoxMax - bBoxMin).two_norm();
186 DUNE_THROW(Dune::InvalidStateException,
"Points do not lie within a plane");
188 auto corners = grahamConvexHull<2>(points);
189 if (corners.size() != 4)
190 DUNE_THROW(Dune::InvalidStateException,
"Points do not span a strictly convex polygon!");
193 std::array<int, 4> orientations;
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.");
199 const auto quadrilateral = GeometryType(Dune::GeometryTypes::quadrilateral, corners);
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");
205 return quadrilateral;
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