23#ifndef DUMUX_GEOMETRY_INTERSECTION_HH
24#define DUMUX_GEOMETRY_INTERSECTION_HH
26#include <dune/common/exceptions.hh>
27#include <dune/common/promotiontraits.hh>
28#include <dune/geometry/referenceelements.hh>
29#include <dune/geometry/multilineargeometry.hh>
37namespace IntersectionPolicy {
40template<
class ct,
int dw>
44 using Point = Dune::FieldVector<ctype, dw>;
52template<
class ct,
int dw>
56 using Point = Dune::FieldVector<ctype, dw>;
65template<
class ct,
int dw>
69 using Point = Dune::FieldVector<ctype, dw>;
78template<
class Geometry1,
class Geometry2,
int dim2 = Geometry2::mydimension>
81 static constexpr int dimworld = Geometry1::coorddimension;
82 static constexpr int isDim = std::min(
int(Geometry1::mydimension),
int(Geometry2::mydimension) );
83 static_assert(dimworld == int(Geometry2::coorddimension),
"Geometries must have the same coordinate dimension!");
84 static_assert(isDim == 1 || isDim == 2,
"No chooser class specialization available for given dimensions");
85 using ctype =
typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
87 using type = std::conditional_t< isDim == 2,
93template<
class Geometry1,
class Geometry2>
113template<
class Geo1,
class Geo2,
class ctype,
114 class GetFacetCornerIndices,
class ComputeNormalFunction >
115bool computeSegmentIntersection(
const Geo1& geo1,
const Geo2& geo2, ctype baseEps,
116 ctype& tfirst, ctype& tlast,
117 const GetFacetCornerIndices& getFacetCornerIndices,
118 const ComputeNormalFunction& computeNormal)
120 static_assert(int(Geo2::mydimension) == 1,
"Geometry2 must be a segment");
121 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
122 "Dimension of Geometry1 must be higher than that of Geometry2");
124 const auto a = geo2.corner(0);
125 const auto b = geo2.corner(1);
126 const auto d = b - a;
134 const auto facets = getFacetCornerIndices(geo1);
135 for (
const auto& f : facets)
137 const auto n = computeNormal(f);
139 const auto c0 = geo1.corner(f[0]);
140 const ctype denom = n*d;
141 const ctype dist = n*(a-c0);
144 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
145 const ctype
eps = baseEps*edge1.two_norm();
150 if (abs(denom) <
eps)
157 const ctype t = -dist / denom;
173 if (tfirst > tlast - baseEps)
191<
class Geometry1,
class Geometry2,
192 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
193 int dimworld = Geometry1::coorddimension,
194 int dim1 = Geometry1::mydimension,
195 int dim2 = Geometry2::mydimension>
198 static constexpr int dimWorld = Policy::dimWorld;
201 using ctype =
typename Policy::ctype;
202 using Point =
typename Policy::Point;
211 static_assert(dimworld == Geometry2::coorddimension,
"Can only intersect geometries of same coordinate dimension");
212 DUNE_THROW(Dune::NotImplemented,
"Geometry intersection detection with intersection computation not implemented for dimworld = "
213 << dimworld <<
", dim1 = " << dim1 <<
", dim2 = " << dim2);
221template <
class Geometry1,
class Geometry2,
class Policy>
224 enum { dimworld = 2 };
229 static constexpr typename Policy::ctype eps_ = 1.5e-7;
232 using ctype =
typename Policy::ctype;
233 using Point =
typename Policy::Point;
245 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
248 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
250 const auto v1 = geo1.corner(1) - geo1.corner(0);
251 const auto v2 = geo2.corner(1) - geo2.corner(0);
252 const auto ac = geo2.corner(0) - geo1.corner(0);
254 auto n2 =
Point({-1.0*v2[1], v2[0]});
258 const auto dist12 = n2*ac;
264 const auto v1Norm2 = v1.two_norm2();
265 const auto eps = eps_*sqrt(v1Norm2);
266 const auto eps2 = eps_*v1Norm2;
268 const auto sp = n2*v1;
271 if (abs(dist12) >
eps)
277 if ( ac.two_norm2() < eps2 )
280 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
285 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
288 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
297 const auto t1 = dist12 / sp;
300 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
304 auto isPoint = geo1.global(t1);
307 const auto c = geo2.corner(0);
308 const auto d = geo2.corner(1);
310 using std::min;
using std::max;
311 std::array<ctype, 4> bBox({ min(c[0], d[0]), min(c[1], d[1]), max(c[0], d[0]), max(c[1], d[1]) });
329 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
332 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
333 DUNE_THROW(Dune::NotImplemented,
"segment-segment intersection detection for segment-like intersections");
341template <
class Geometry1,
class Geometry2,
class Policy>
344 enum { dimworld = 2 };
349 using ctype =
typename Policy::ctype;
350 using Point =
typename Policy::Point;
357 static constexpr ctype eps_ = 1.5e-7;
358 using ReferenceElements =
typename Dune::ReferenceElements<ctype, dim1>;
368 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
371 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
374 if (intersect_(geo1, geo2, tfirst, tlast))
392 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
395 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
398 if (!intersect_(geo1, geo2, tfirst, tlast))
402 if ( tfirst > tlast - eps_ )
416 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
419 auto getFacetCorners = [] (
const Geometry1& geo1)
421 std::vector< std::array<int, 2> > facetCorners;
422 switch (geo1.corners())
425 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
428 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
431 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
432 << geo1.type() <<
" with "<< geo1.corners() <<
" corners.");
439 const auto center1 = geo1.center();
440 auto computeNormal = [&geo1, ¢er1] (
const std::array<int, 2>& facetCorners)
442 const auto c0 = geo1.corner(facetCorners[0]);
443 const auto c1 = geo1.corner(facetCorners[1]);
444 const auto edge = c1 - c0;
446 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
451 if ( n*(center1-c0) > 0.0 )
457 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
465template <
class Geometry1,
class Geometry2,
class Policy>
479 template<
class P = Policy>
490template <
class Geometry1,
class Geometry2,
class Policy>
493 enum { dimworld = 2 };
498 using ctype =
typename Policy::ctype;
499 using Point =
typename Policy::Point;
506 static constexpr ctype eps_ = 1.5e-7;
507 using ReferenceElementsGeo1 =
typename Dune::ReferenceElements<ctype, dim1>;
508 using ReferenceElementsGeo2 =
typename Dune::ReferenceElements<ctype, dim2>;
523 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
526 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
529 std::vector<Point> points; points.reserve(6);
532 for (
int i = 0; i < geo1.corners(); ++i)
534 points.emplace_back(geo1.corner(i));
536 const auto numPoints1 = points.size();
537 if (numPoints1 != geo1.corners())
540 for (
int i = 0; i < geo2.corners(); ++i)
542 points.emplace_back(geo2.corner(i));
547 if (points.size() - numPoints1 != geo2.corners())
549 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
550 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
553 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
555 for (
int i = 0; i < referenceElement1.size(dim1-1); ++i)
557 const auto localEdgeGeom1 = referenceElement1.template geometry<dim1-1>(i);
558 const auto edge1 = SegGeometry( Dune::GeometryTypes::line,
559 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
560 geo1.global(localEdgeGeom1.corner(1))} ));
562 for (
int j = 0; j < referenceElement2.size(dim2-1); ++j)
564 const auto localEdgeGeom2 = referenceElement2.template geometry<dim2-1>(j);
565 const auto edge2 = SegGeometry( Dune::GeometryTypes::line,
566 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
567 geo2.global(localEdgeGeom2.corner(1))} ));
570 typename EdgeTest::Intersection edgeIntersection;
571 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
572 points.emplace_back(edgeIntersection);
582 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
583 std::sort(points.begin(), points.end(), [&
eps](
const auto& a,
const auto& b) ->
bool
586 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
589 auto removeIt = std::unique(points.begin(), points.end(), [&
eps](
const auto& a,
const auto&b)
591 return (b-a).two_norm() < eps;
594 points.erase(removeIt, points.end());
597 if (points.size() < 3)
613 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
616 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
617 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for segment-like intersections");
627 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
630 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
631 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for touching points");
639template <
class Geometry1,
class Geometry2,
class Policy>
642 enum { dimworld = 3 };
647 using ctype =
typename Policy::ctype;
648 using Point =
typename Policy::Point;
655 static constexpr ctype eps_ = 1.5e-7;
656 using ReferenceElements =
typename Dune::ReferenceElements<ctype, dim1>;
670 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
673 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
676 if (intersect_(geo1, geo2, tfirst, tlast))
698 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
701 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
704 if (intersect_(geo1, geo2, tfirst, tlast))
708 if ( abs(tfirst - tlast) < eps_ )
722 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
725 auto getFacetCorners = [] (
const Geometry1& geo1)
727 std::vector< std::vector<int> > facetCorners;
729 switch (geo1.corners())
733 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
734 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
737 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
740 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
741 << geo1.type() <<
", "<< geo1.corners() <<
" corners.");
748 auto computeNormal = [&geo1] (
const std::vector<int>& facetCorners)
750 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
751 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
759 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
767template <
class Geometry1,
class Geometry2,
class Policy>
780 template<
class P = Policy>
791template <
class Geometry1,
class Geometry2,
class Policy>
794 enum { dimworld = 3 };
799 using ctype =
typename Policy::ctype;
800 using Point =
typename Policy::Point;
807 static constexpr ctype eps_ = 1.5e-7;
808 using ReferenceElementsGeo1 =
typename Dune::ReferenceElements<ctype, dim1>;
809 using ReferenceElementsGeo2 =
typename Dune::ReferenceElements<ctype, dim2>;
826 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
829 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
832 std::vector<Point> points; points.reserve(10);
835 for (
int i = 0; i < geo1.corners(); ++i)
837 points.emplace_back(geo1.corner(i));
840 for (
int i = 0; i < geo2.corners(); ++i)
842 points.emplace_back(geo2.corner(i));
845 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
846 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
848 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
849 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
855 for (
int i = 0; i < referenceElement1.size(dim1-1); ++i)
857 const auto localEdgeGeom = referenceElement1.template geometry<dim1-1>(i);
858 const auto p = geo1.global(localEdgeGeom.corner(0));
859 const auto q = geo1.global(localEdgeGeom.corner(1));
860 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
863 typename PolySegTest::Intersection polySegIntersection;
864 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
865 points.emplace_back(polySegIntersection);
869 for (
int i = 0; i < referenceElement1.size(1); ++i)
871 const auto faceGeo = [&]()
873 const auto localFaceGeo = referenceElement1.template geometry<1>(i);
874 if (localFaceGeo.corners() == 4)
876 const auto a = geo1.global(localFaceGeo.corner(0));
877 const auto b = geo1.global(localFaceGeo.corner(1));
878 const auto c = geo1.global(localFaceGeo.corner(2));
879 const auto d = geo1.global(localFaceGeo.corner(3));
881 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
885 const auto a = geo1.global(localFaceGeo.corner(0));
886 const auto b = geo1.global(localFaceGeo.corner(1));
887 const auto c = geo1.global(localFaceGeo.corner(2));
889 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
893 for (
int j = 0; j < referenceElement2.size(1); ++j)
895 const auto localEdgeGeom = referenceElement2.template geometry<1>(j);
896 const auto p = geo2.global(localEdgeGeom.corner(0));
897 const auto q = geo2.global(localEdgeGeom.corner(1));
899 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
902 typename PolySegTest::Intersection polySegIntersection;
903 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
904 points.emplace_back(polySegIntersection);
909 if (points.empty())
return false;
912 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
913 std::sort(points.begin(), points.end(), [&
eps](
const auto& a,
const auto& b) ->
bool
916 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])));
919 auto removeIt = std::unique(points.begin(), points.end(), [&
eps](
const auto& a,
const auto&b)
921 return (b-a).two_norm() < eps;
924 points.erase(removeIt, points.end());
927 if (points.size() < 3)
return false;
950 template<
class P = Policy, std::enable_if_t<P::dimIntersection != 2,
int> = 0>
953 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
954 DUNE_THROW(Dune::NotImplemented,
"Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
962template <
class Geometry1,
class Geometry2,
class Policy>
975 template<
class P = Policy>
986template <
class Geometry1,
class Geometry2,
class Policy>
989 enum { dimworld = 3 };
994 using ctype =
typename Policy::ctype;
995 using Point =
typename Policy::Point;
1002 static constexpr ctype eps_ = 1.5e-7;
1003 using ReferenceElements =
typename Dune::ReferenceElements<ctype, dim1>;
1015 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1018 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1020 ctype tfirst, tlast;
1021 if (intersect_(geo1, geo2, tfirst, tlast))
1041 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1044 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1046 const auto p = geo2.corner(0);
1047 const auto q = geo2.corner(1);
1049 const auto a = geo1.corner(0);
1050 const auto b = geo1.corner(1);
1051 const auto c = geo1.corner(2);
1053 if (geo1.corners() == 3)
1054 return intersect_<Policy>(a, b, c, p, q, is);
1056 else if (geo1.corners() == 4)
1059 bool hasSegment1, hasSegment2;
1061 const auto d = geo1.corner(3);
1062 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1063 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1065 if (hasSegment1 || hasSegment2)
1068 if (intersects1) { is = is1;
return true; }
1069 if (intersects2) { is = is2;
return true; }
1075 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
1076 << geo1.type() <<
", "<< geo1.corners() <<
" corners.");
1089 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1090 [[deprecated(
"Please use intersection(triangle, segment, ...) instead")]]
1095 return intersect_<Policy>(a, b, c, p, q, is);
1102 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1103 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
1106 auto getFacetCorners = [] (
const Geometry1& geo1)
1108 std::vector< std::array<int, 2> > facetCorners;
1109 switch (geo1.corners())
1112 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1115 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1118 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
1119 << geo1.type() <<
" with "<< geo1.corners() <<
" corners.");
1122 return facetCorners;
1125 const auto center1 = geo1.center();
1126 const auto normal1 =
crossProduct(geo1.corner(1) - geo1.corner(0),
1127 geo1.corner(2) - geo1.corner(0));
1130 auto computeNormal = [¢er1, &normal1, &geo1] (
const std::array<int, 2>& facetCorners)
1132 const auto c0 = geo1.corner(facetCorners[0]);
1133 const auto c1 = geo1.corner(facetCorners[1]);
1134 const auto edge = c1 - c0;
1136 Dune::FieldVector<ctype, dimworld> n =
crossProduct(edge, normal1);
1141 if ( n*(center1-c0) > 0.0 )
1147 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
1153 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1154 static bool intersect_(
const Point& a,
const Point& b,
const Point& c,
1159 return intersect_(a, b, c, p, q, is, hasSegment);
1166 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1167 static bool intersect_(
const Point& a,
const Point& b,
const Point& c,
1171 hasSegmentIs =
false;
1173 const auto ab = b - a;
1174 const auto ac = c - a;
1175 const auto qp = p - q;
1182 const auto denom = normal*qp;
1183 const auto abnorm2 = ab.two_norm2();
1184 const auto eps = eps_*abnorm2*qp.two_norm();
1186 if (abs(denom) <
eps)
1188 const auto pa = a - p;
1189 const auto denom2 = normal*pa;
1190 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1195 using SegmentPolicy =
typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1196 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1197 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1198 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1199 using SegmentIntersectionType =
typename SegmentIntersectionAlgorithm::Intersection;
1200 SegmentIntersectionType segmentIs;
1202 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1203 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1204 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1206 hasSegmentIs =
true;
1212 for (
const auto& ip : {p, q})
1228 const auto ap = p - a;
1229 const auto t = (ap*normal)/denom;
1230 if (t < 0.0 - eps_)
return false;
1231 if (t > 1.0 + eps_)
return false;
1236 const auto v = (ac*e)/denom;
1237 if (v < -eps_ || v > 1.0 + eps_)
return false;
1238 const auto w = -(ab*e)/denom;
1239 if (w < -eps_ || v + w > 1.0 + eps_)
return false;
1243 const auto u = 1.0 - v - w;
1258template <
class Geometry1,
class Geometry2,
class Policy>
1271 template<
class P = Policy>
A function to compute the convex hull of a point cloud and a function to triangulate the polygon area...
Detect if a point intersects a geometry.
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
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:38
bool intersectsPointSimplex(const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &p0, const Dune::FieldVector< ctype, dimworld > &p1, const Dune::FieldVector< ctype, dimworld > &p2, const Dune::FieldVector< ctype, dimworld > &p3)
Find out whether a point is inside a tetrahedron (p0, p1, p2, p3) (dimworld is 3)
Definition: intersectspointsimplex.hh:36
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
bool intersectsPointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Check whether a point is intersectin a bounding box (dimworld == 3)
Definition: geometry/boundingboxtree.hh:304
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:94
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
Policy structure for point-like intersections.
Definition: geometryintersection.hh:42
static constexpr int dimIntersection
Definition: geometryintersection.hh:47
Point Intersection
Definition: geometryintersection.hh:48
ct ctype
Definition: geometryintersection.hh:43
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:44
static constexpr int dimWorld
Definition: geometryintersection.hh:46
Policy structure for segment-like intersections.
Definition: geometryintersection.hh:54
Segment Intersection
Definition: geometryintersection.hh:61
static constexpr int dimIntersection
Definition: geometryintersection.hh:60
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:56
static constexpr int dimWorld
Definition: geometryintersection.hh:59
std::array< Point, 2 > Segment
Definition: geometryintersection.hh:57
ct ctype
Definition: geometryintersection.hh:55
Policy structure for polygon-like intersections.
Definition: geometryintersection.hh:67
static constexpr int dimWorld
Definition: geometryintersection.hh:72
static constexpr int dimIntersection
Definition: geometryintersection.hh:73
ct ctype
Definition: geometryintersection.hh:68
Polygon Intersection
Definition: geometryintersection.hh:74
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:69
std::vector< Point > Polygon
Definition: geometryintersection.hh:70
default policy chooser class
Definition: geometryintersection.hh:80
std::conditional_t< isDim==2, PolygonPolicy< ctype, Geometry1::coorddimension >, SegmentPolicy< ctype, Geometry1::coorddimension > > type
Definition: geometryintersection.hh:89
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:197
typename Policy::ctype ctype
Definition: geometryintersection.hh:201
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:203
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:206
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Determine if the two geometries intersect and compute the intersection geometry.
Definition: geometryintersection.hh:209
typename Policy::Point Point
Definition: geometryintersection.hh:202
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:237
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:234
typename Policy::ctype ctype
Definition: geometryintersection.hh:232
typename Policy::Point Point
Definition: geometryintersection.hh:233
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:246
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:369
typename Policy::ctype ctype
Definition: geometryintersection.hh:349
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:354
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:351
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:393
typename Policy::Point Point
Definition: geometryintersection.hh:350
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:480
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:500
typename Policy::ctype ctype
Definition: geometryintersection.hh:498
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:503
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:524
typename Policy::Point Point
Definition: geometryintersection.hh:499
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:649
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:652
typename Policy::Point Point
Definition: geometryintersection.hh:648
typename Policy::ctype ctype
Definition: geometryintersection.hh:647
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:671
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:781
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:801
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:804
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:827
typename Policy::Point Point
Definition: geometryintersection.hh:800
typename Policy::ctype ctype
Definition: geometryintersection.hh:799
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:976
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:999
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:996
typename Policy::Point Point
Definition: geometryintersection.hh:995
static bool intersection(const Point &a, const Point &b, const Point &c, const Point &p, const Point &q, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1091
typename Policy::ctype ctype
Definition: geometryintersection.hh:994
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1016
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1042
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:1272
An axis-aligned bounding box volume hierarchy for dune grids.