12#ifndef DUMUX_GEOMETRY_GRAHAM_CONVEX_HULL_HH 
   13#define DUMUX_GEOMETRY_GRAHAM_CONVEX_HULL_HH 
   20#include <dune/common/exceptions.hh> 
   21#include <dune/common/fvector.hh> 
   36                   const Dune::FieldVector<ctype, 3>& b,
 
   37                   const Dune::FieldVector<ctype, 3>& c,
 
   38                   const Dune::FieldVector<ctype, 3>& 
normal)
 
   43    const auto area = f*
normal;
 
 
   54template<
int dim, 
class ctype,
 
   55         std::enable_if_t<(dim==2), 
int> = 0>
 
   56std::vector<Dune::FieldVector<ctype, 3>>
 
   59    using Point = Dune::FieldVector<ctype, 3>;
 
   60    std::vector<Point> convexHull;
 
   63    if (points.size() < 3)
 
   67    if (points.size() == 3)
 
   71    const auto a = points[1] - points[0];
 
   72    auto b = points[2] - points[0];
 
   77    auto norm = 
normal.two_norm();
 
   78    while (norm == 0.0 && k < points.size()-1)
 
   90    const auto eps = 1e-7*sqrt(norm);
 
   94    auto minIt = std::min_element(points.begin(), points.end(), [&eps](
const auto& a, 
const auto& b)
 
   97        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])));
 
  101    std::iter_swap(minIt, points.begin());
 
  105    const auto pivot = points[0];
 
  106    std::sort(points.begin()+1, points.end(), [&](
const auto& a, 
const auto& b)
 
  108        const int order = getOrientation(pivot, a, b, normal);
 
  110             return (a-pivot).two_norm() < (b-pivot).two_norm();
 
  112             return (order == -1);
 
  116    convexHull.reserve(50);
 
  117    convexHull.push_back(points[0]);
 
  118    convexHull.push_back(points[1]);
 
  119    convexHull.push_back(points[2]);
 
  124    for (std::size_t i = 3; i < points.size(); ++i)
 
  126        Point p = convexHull.back();
 
  127        convexHull.pop_back();
 
  132            if (convexHull.size() == 1)
 
  136                assert(i < points.size()-1);
 
  141                p = convexHull.back();
 
  142                convexHull.pop_back();
 
  147        convexHull.emplace_back(std::move(p));
 
  148        convexHull.push_back(points[i]);
 
 
  160template<
int dim, 
class ctype,
 
  161         std::enable_if_t<(dim==2), 
int> = 0>
 
  162std::vector<Dune::FieldVector<ctype, 2>>
 
  165    std::vector<Dune::FieldVector<ctype, 3>> points3D;
 
  166    points3D.reserve(points.size());
 
  167    std::transform(points.begin(), points.end(), std::back_inserter(points3D),
 
  168                   [](
const auto& p) { return Dune::FieldVector<ctype, 3>({p[0], p[1], 0.0}); });
 
  172    std::vector<Dune::FieldVector<ctype, 2>> result2D;
 
  173    result2D.reserve(result3D.size());
 
  174    std::transform(result3D.begin(), result3D.end(), std::back_inserter(result2D),
 
  175                   [](
const auto& p) { return Dune::FieldVector<ctype, 2>({p[0], p[1]}); });
 
 
  185template<
int dim, 
class ctype, 
int dimWorld>
 
  186std::vector<Dune::FieldVector<ctype, dimWorld>> 
grahamConvexHull(std::vector<Dune::FieldVector<ctype, dimWorld>>& points)
 
 
  197template<
int dim, 
class ctype, 
int dimWorld>
 
  198std::vector<Dune::FieldVector<ctype, dimWorld>> 
grahamConvexHull(
const std::vector<Dune::FieldVector<ctype, dimWorld>>& points)
 
  200    auto copyPoints = points;
 
 
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
constexpr int sign(const ValueType &value) noexcept
Sign or signum function.
Definition math.hh:658
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:57
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< 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.