12#ifndef DUMUX_GEOMETRY_MAKE_GEOMETRY_HH 
   13#define DUMUX_GEOMETRY_MAKE_GEOMETRY_HH 
   18#include <dune/common/fvector.hh> 
   19#include <dune/common/fmatrix.hh> 
   20#include <dune/common/exceptions.hh> 
   21#include <dune/geometry/multilineargeometry.hh> 
   31template<
class CoordScalar>
 
   32bool pointsAreCoplanar(
const std::vector<Dune::FieldVector<CoordScalar, 3>>& points, 
const CoordScalar scale)
 
   34    if (points.size() != 4)
 
   35        DUNE_THROW(Dune::InvalidStateException, 
"Check only works for 4 points!");
 
   38    Dune::FieldMatrix<CoordScalar, 4, 4> M;
 
   39    for (
int i = 0; i < 3; ++i )
 
   40        M[i] = {points[0][i], points[1][i], points[2][i], points[3][i]};
 
   41    M[3] = {1.0*scale, 1.0*scale, 1.0*scale, 1.0*scale};
 
   44    return abs(M.determinant()) < 1.5e-7*scale*scale*scale*scale;
 
 
   51template<
class CoordScalar>
 
   54    Dune::FieldVector<CoordScalar, 3> bBoxMin(std::numeric_limits<CoordScalar>::max());
 
   55    Dune::FieldVector<CoordScalar, 3> bBoxMax(std::numeric_limits<CoordScalar>::lowest());
 
   56    for (
const auto& p : points)
 
   58        for (
int i=0; i<3; i++)
 
   62            bBoxMin[i] = min(bBoxMin[i], p[i]);
 
   63            bBoxMax[i] = max(bBoxMax[i], p[i]);
 
   67    const auto size = (bBoxMax - bBoxMin).two_norm();
 
 
   79template<
class CoordScalar>
 
   80std::vector<Dune::FieldVector<CoordScalar, 3>> 
getReorderedPoints(
const std::vector<Dune::FieldVector<CoordScalar, 3>>& points)
 
   82    std::array<int, 4> tmp;
 
 
   93template<
class CoordScalar>
 
   94std::vector<Dune::FieldVector<CoordScalar, 3>> 
getReorderedPoints(
const std::vector<Dune::FieldVector<CoordScalar, 3>>& points,
 
   95                                                             std::array<int, 4>& orientations)
 
   97    if(points.size() == 4)
 
  100        auto& p1 = points[1];
 
  101        auto& p2 = points[2];
 
  102        auto& p3 = points[3];
 
  114        const bool diagonalsIntersect = (orientations[0] != orientations[1]) && (orientations[2] != orientations[3]);
 
  117        if(diagonalsIntersect)
 
  121        using GlobalPosition = Dune::FieldVector<CoordScalar, 3>;
 
  122        if(orientations[0] == 1)
 
  123            return std::vector<GlobalPosition>{p1, p0, p2, p3};
 
  124        else if(orientations[0] == -1)
 
  125            return std::vector<GlobalPosition>{p3, p1, p0, p2};
 
  127            DUNE_THROW(Dune::InvalidStateException, 
"Could not reorder points");
 
  130        DUNE_THROW(Dune::NotImplemented, 
"Reorder for " << points.size() << 
" points.");
 
 
  141template<
class CoordScalar, 
bool enableSanityCheck = true>
 
  144    if (points.size() != 4)
 
  145        DUNE_THROW(Dune::InvalidStateException, 
"A quadrilateral needs 4 corner points!");
 
  147    using GlobalPosition = Dune::FieldVector<CoordScalar, 3>;
 
  148    static constexpr auto coordDim = GlobalPosition::dimension;
 
  149    static constexpr auto dim = coordDim-1;
 
  150    using GeometryType = Dune::MultiLinearGeometry<CoordScalar, dim, coordDim>;
 
  153    if (!enableSanityCheck)
 
  154        return GeometryType(Dune::GeometryTypes::quadrilateral, points);
 
  157    Dune::FieldVector<CoordScalar, 3> bBoxMin(std::numeric_limits<CoordScalar>::max());
 
  158    Dune::FieldVector<CoordScalar, 3> bBoxMax(std::numeric_limits<CoordScalar>::lowest());
 
  159    for (
const auto& p : points)
 
  161        for (
int i = 0; i < 3; i++)
 
  165            bBoxMin[i] = min(bBoxMin[i], p[i]);
 
  166            bBoxMax[i] = max(bBoxMax[i], p[i]);
 
  170    const auto size = (bBoxMax - bBoxMin).two_norm();
 
  174        DUNE_THROW(Dune::InvalidStateException, 
"Points do not lie within a plane");
 
  177    if (corners.size() != 4)
 
  178        DUNE_THROW(Dune::InvalidStateException, 
"Points do not span a strictly convex polygon!");
 
  181    std::array<int, 4> orientations;
 
  184    if (std::any_of(orientations.begin(), orientations.end(), [](
auto i){ return i == 0; }))
 
  185        DUNE_THROW(Dune::InvalidStateException, 
"More than two points lie on the same line.");
 
  187    const auto quadrilateral = GeometryType(Dune::GeometryTypes::quadrilateral, corners);
 
  189    const auto eps = 1e-7;
 
  190    if (quadrilateral.volume() < eps*size*size)
 
  191        DUNE_THROW(Dune::InvalidStateException, 
"Something went wrong, geometry has area of zero");
 
  193    return quadrilateral;
 
 
A function to compute the convex hull of a point cloud.
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:671
auto makeDuneQuadrilaterial(const std::vector< Dune::FieldVector< CoordScalar, 3 > > &points)
Creates a dune quadrilateral geometry given 4 corner points.
Definition makegeometry.hh:142
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:35
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition normal.hh:26
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:80
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:32
std::vector< Dune::FieldVector< ctype, dimWorld > > grahamConvexHull(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:186
Define some often used mathematical functions.