23#ifndef DUMUX_GEOMETRY_INTERSECTION_HH
24#define DUMUX_GEOMETRY_INTERSECTION_HH
28#include <dune/common/exceptions.hh>
29#include <dune/common/promotiontraits.hh>
30#include <dune/geometry/referenceelements.hh>
31#include <dune/geometry/multilineargeometry.hh>
39namespace IntersectionPolicy {
42template<
class ct,
int dw>
46 using Point = Dune::FieldVector<ctype, dw>;
54template<
class ct,
int dw>
58 using Point = Dune::FieldVector<ctype, dw>;
67template<
class ct,
int dw>
71 using Point = Dune::FieldVector<ctype, dw>;
80template<
class ct,
int dw>
84 using Point = Dune::FieldVector<ctype, dw>;
93template<
class Geometry1,
class Geometry2>
96 static constexpr int dimworld = Geometry1::coorddimension;
97 static constexpr int isDim = std::min(
int(Geometry1::mydimension),
int(Geometry2::mydimension) );
98 static_assert(dimworld == int(Geometry2::coorddimension),
99 "Geometries must have the same coordinate dimension!");
100 static_assert(int(Geometry1::mydimension) <= 3 && int(Geometry2::mydimension) <= 3,
101 "Geometries must have dimension 3 or less.");
102 using ctype =
typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
104 using DefaultPolicies = std::tuple<PointPolicy<ctype, dimworld>,
109 using type = std::tuple_element_t<isDim, DefaultPolicies>;
113template<
class Geometry1,
class Geometry2>
133template<
class Geo1,
class Geo2,
class ctype,
134 class GetFacetCornerIndices,
class ComputeNormalFunction >
135bool computeSegmentIntersection(
const Geo1& geo1,
const Geo2& geo2, ctype baseEps,
136 ctype& tfirst, ctype& tlast,
137 const GetFacetCornerIndices& getFacetCornerIndices,
138 const ComputeNormalFunction& computeNormal)
140 static_assert(int(Geo2::mydimension) == 1,
"Geometry2 must be a segment");
141 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
142 "Dimension of Geometry1 must be higher than that of Geometry2");
144 const auto a = geo2.corner(0);
145 const auto b = geo2.corner(1);
146 const auto d = b - a;
154 const auto facets = getFacetCornerIndices(geo1);
155 for (
const auto& f : facets)
157 const auto n = computeNormal(f);
159 const auto c0 = geo1.corner(f[0]);
160 const ctype denom = n*d;
161 const ctype dist = n*(a-c0);
164 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
165 const ctype eps = baseEps*edge1.two_norm();
170 if (abs(denom) < eps)
177 const ctype t = -dist / denom;
193 if (tfirst > tlast - baseEps)
211<
class Geometry1,
class Geometry2,
212 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
213 int dimworld = Geometry1::coorddimension,
214 int dim1 = Geometry1::mydimension,
215 int dim2 = Geometry2::mydimension>
218 static constexpr int dimWorld = Policy::dimWorld;
221 using ctype =
typename Policy::ctype;
222 using Point =
typename Policy::Point;
228 static_assert(dimworld == Geometry2::coorddimension,
"Can only intersect geometries of same coordinate dimension");
229 DUNE_THROW(Dune::NotImplemented,
"Geometry intersection detection with intersection computation not implemented for dimworld = "
230 << dimworld <<
", dim1 = " << dim1 <<
", dim2 = " << dim2);
238template <
class Geometry1,
class Geometry2,
class Policy>
241 enum { dimworld = 2 };
246 static constexpr typename Policy::ctype eps_ = 1.5e-7;
249 using ctype =
typename Policy::ctype;
250 using Point =
typename Policy::Point;
259 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
262 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
264 const auto v1 = geo1.corner(1) - geo1.corner(0);
265 const auto v2 = geo2.corner(1) - geo2.corner(0);
266 const auto ac = geo2.corner(0) - geo1.corner(0);
268 auto n2 =
Point({-1.0*v2[1], v2[0]});
272 const auto dist12 = n2*ac;
278 const auto v1Norm2 = v1.two_norm2();
279 const auto eps = eps_*sqrt(v1Norm2);
280 const auto eps2 = eps_*v1Norm2;
282 const auto sp = n2*v1;
285 if (abs(dist12) > eps)
291 if ( ac.two_norm2() < eps2 )
294 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
299 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
302 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
311 const auto t1 = dist12 / sp;
314 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
318 auto isPoint = geo1.global(t1);
321 const auto c = geo2.corner(0);
322 const auto d = geo2.corner(1);
324 using std::min;
using std::max;
325 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]) });
343 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
346 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
347 DUNE_THROW(Dune::NotImplemented,
"segment-segment intersection detection for segment-like intersections");
355template <
class Geometry1,
class Geometry2,
class Policy>
358 enum { dimworld = 2 };
363 using ctype =
typename Policy::ctype;
364 using Point =
typename Policy::Point;
368 static constexpr ctype eps_ = 1.5e-7;
369 using ReferenceElements =
typename Dune::ReferenceElements<ctype, dim1>;
379 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
382 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
385 if (intersect_(geo1, geo2, tfirst, tlast))
403 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
406 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
409 if (!intersect_(geo1, geo2, tfirst, tlast))
413 if ( tfirst > tlast - eps_ )
427 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
430 auto getFacetCorners = [] (
const Geometry1& geo1)
432 std::vector< std::array<int, 2> > facetCorners;
433 switch (geo1.corners())
436 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
439 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
442 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
443 << geo1.type() <<
" with "<< geo1.corners() <<
" corners.");
450 const auto center1 = geo1.center();
451 auto computeNormal = [&geo1, ¢er1] (
const std::array<int, 2>& facetCorners)
453 const auto c0 = geo1.corner(facetCorners[0]);
454 const auto c1 = geo1.corner(facetCorners[1]);
455 const auto edge = c1 - c0;
457 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
462 if ( n*(center1-c0) > 0.0 )
468 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
476template <
class Geometry1,
class Geometry2,
class Policy>
490 template<
class P = Policy>
501template <
class Geometry1,
class Geometry2,
class Policy>
504 enum { dimworld = 2 };
509 using ctype =
typename Policy::ctype;
510 using Point =
typename Policy::Point;
514 static constexpr ctype eps_ = 1.5e-7;
515 using ReferenceElementsGeo1 =
typename Dune::ReferenceElements<ctype, dim1>;
516 using ReferenceElementsGeo2 =
typename Dune::ReferenceElements<ctype, dim2>;
531 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
534 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
537 std::vector<Point> points; points.reserve(6);
540 for (
int i = 0; i < geo1.corners(); ++i)
542 points.emplace_back(geo1.corner(i));
544 const auto numPoints1 = points.size();
545 if (numPoints1 != geo1.corners())
548 for (
int i = 0; i < geo2.corners(); ++i)
550 points.emplace_back(geo2.corner(i));
555 if (points.size() - numPoints1 != geo2.corners())
557 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
558 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
561 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
563 for (
int i = 0; i < referenceElement1.size(dim1-1); ++i)
565 const auto localEdgeGeom1 = referenceElement1.template geometry<dim1-1>(i);
566 const auto edge1 = SegGeometry( Dune::GeometryTypes::line,
567 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
568 geo1.global(localEdgeGeom1.corner(1))} ));
570 for (
int j = 0; j < referenceElement2.size(dim2-1); ++j)
572 const auto localEdgeGeom2 = referenceElement2.template geometry<dim2-1>(j);
573 const auto edge2 = SegGeometry( Dune::GeometryTypes::line,
574 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
575 geo2.global(localEdgeGeom2.corner(1))} ));
578 typename EdgeTest::Intersection edgeIntersection;
579 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
580 points.emplace_back(edgeIntersection);
590 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
591 std::sort(points.begin(), points.end(), [&eps](
const auto& a,
const auto& b) ->
bool
594 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
597 auto removeIt = std::unique(points.begin(), points.end(), [&eps](
const auto& a,
const auto&b)
599 return (b-a).two_norm() < eps;
602 points.erase(removeIt, points.end());
605 if (points.size() < 3)
621 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
624 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
625 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for segment-like intersections");
635 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
638 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
639 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for touching points");
647template <
class Geometry1,
class Geometry2,
class Policy>
650 enum { dimworld = 3 };
655 using ctype =
typename Policy::ctype;
656 using Point =
typename Policy::Point;
660 static constexpr ctype eps_ = 1.5e-7;
661 using ReferenceElements =
typename Dune::ReferenceElements<ctype, dim1>;
675 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
678 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
681 if (intersect_(geo1, geo2, tfirst, tlast))
703 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
706 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
709 if (intersect_(geo1, geo2, tfirst, tlast))
713 if ( abs(tfirst - tlast) < eps_ )
727 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
730 auto getFacetCorners = [] (
const Geometry1& geo1)
732 std::vector< std::vector<int> > facetCorners;
734 switch (geo1.corners())
738 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
739 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
742 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
745 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
746 << geo1.type() <<
", "<< geo1.corners() <<
" corners.");
753 auto computeNormal = [&geo1] (
const std::vector<int>& facetCorners)
755 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
756 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
764 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
772template <
class Geometry1,
class Geometry2,
class Policy>
785 template<
class P = Policy>
796template <
class Geometry1,
class Geometry2,
class Policy>
799 enum { dimworld = 3 };
804 using ctype =
typename Policy::ctype;
805 using Point =
typename Policy::Point;
809 static constexpr ctype eps_ = 1.5e-7;
810 using ReferenceElementsGeo1 =
typename Dune::ReferenceElements<ctype, dim1>;
811 using ReferenceElementsGeo2 =
typename Dune::ReferenceElements<ctype, dim2>;
828 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
831 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
834 std::vector<Point> points; points.reserve(10);
837 for (
int i = 0; i < geo1.corners(); ++i)
839 points.emplace_back(geo1.corner(i));
842 for (
int i = 0; i < geo2.corners(); ++i)
844 points.emplace_back(geo2.corner(i));
847 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
848 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
850 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
851 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
857 for (
int i = 0; i < referenceElement1.size(dim1-1); ++i)
859 const auto localEdgeGeom = referenceElement1.template geometry<dim1-1>(i);
860 const auto p = geo1.global(localEdgeGeom.corner(0));
861 const auto q = geo1.global(localEdgeGeom.corner(1));
862 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
865 typename PolySegTest::Intersection polySegIntersection;
866 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
867 points.emplace_back(polySegIntersection);
871 for (
int i = 0; i < referenceElement1.size(1); ++i)
873 const auto faceGeo = [&]()
875 const auto localFaceGeo = referenceElement1.template geometry<1>(i);
876 if (localFaceGeo.corners() == 4)
878 const auto a = geo1.global(localFaceGeo.corner(0));
879 const auto b = geo1.global(localFaceGeo.corner(1));
880 const auto c = geo1.global(localFaceGeo.corner(2));
881 const auto d = geo1.global(localFaceGeo.corner(3));
883 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
887 const auto a = geo1.global(localFaceGeo.corner(0));
888 const auto b = geo1.global(localFaceGeo.corner(1));
889 const auto c = geo1.global(localFaceGeo.corner(2));
891 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
895 for (
int j = 0; j < referenceElement2.size(1); ++j)
897 const auto localEdgeGeom = referenceElement2.template geometry<1>(j);
898 const auto p = geo2.global(localEdgeGeom.corner(0));
899 const auto q = geo2.global(localEdgeGeom.corner(1));
901 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
904 typename PolySegTest::Intersection polySegIntersection;
905 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
906 points.emplace_back(polySegIntersection);
911 if (points.empty())
return false;
914 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
915 std::sort(points.begin(), points.end(), [&eps](
const auto& a,
const auto& b) ->
bool
918 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])));
921 auto removeIt = std::unique(points.begin(), points.end(), [&eps](
const auto& a,
const auto&b)
923 return (b-a).two_norm() < eps;
926 points.erase(removeIt, points.end());
929 if (points.size() < 3)
return false;
952 template<
class P = Policy, std::enable_if_t<P::dimIntersection != 2,
int> = 0>
955 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
956 DUNE_THROW(Dune::NotImplemented,
"Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
964template <
class Geometry1,
class Geometry2,
class Policy>
977 template<
class P = Policy>
988template <
class Geometry1,
class Geometry2,
class Policy>
991 enum { dimworld = 3 };
996 using ctype =
typename Policy::ctype;
997 using Point =
typename Policy::Point;
1001 static constexpr ctype eps_ = 1.5e-7;
1002 using ReferenceElements =
typename Dune::ReferenceElements<ctype, dim1>;
1014 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1017 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1019 ctype tfirst, tlast;
1020 if (intersect_(geo1, geo2, tfirst, tlast))
1040 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1043 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1045 const auto p = geo2.corner(0);
1046 const auto q = geo2.corner(1);
1048 const auto a = geo1.corner(0);
1049 const auto b = geo1.corner(1);
1050 const auto c = geo1.corner(2);
1052 if (geo1.corners() == 3)
1053 return intersect_<Policy>(a, b, c, p, q, is);
1055 else if (geo1.corners() == 4)
1058 bool hasSegment1, hasSegment2;
1060 const auto d = geo1.corner(3);
1061 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1062 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1064 if (hasSegment1 || hasSegment2)
1067 if (intersects1) { is = is1;
return true; }
1068 if (intersects2) { is = is2;
return true; }
1074 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
1075 << geo1.type() <<
", "<< geo1.corners() <<
" corners.");
1082 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1083 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
1086 auto getFacetCorners = [] (
const Geometry1& geo1)
1088 std::vector< std::array<int, 2> > facetCorners;
1089 switch (geo1.corners())
1092 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1095 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1098 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
1099 << geo1.type() <<
" with "<< geo1.corners() <<
" corners.");
1102 return facetCorners;
1105 const auto center1 = geo1.center();
1106 const auto normal1 =
crossProduct(geo1.corner(1) - geo1.corner(0),
1107 geo1.corner(2) - geo1.corner(0));
1110 auto computeNormal = [¢er1, &normal1, &geo1] (
const std::array<int, 2>& facetCorners)
1112 const auto c0 = geo1.corner(facetCorners[0]);
1113 const auto c1 = geo1.corner(facetCorners[1]);
1114 const auto edge = c1 - c0;
1116 Dune::FieldVector<ctype, dimworld> n =
crossProduct(edge, normal1);
1121 if ( n*(center1-c0) > 0.0 )
1127 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
1133 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1134 static bool intersect_(
const Point& a,
const Point& b,
const Point& c,
1139 return intersect_(a, b, c, p, q, is, hasSegment);
1146 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1147 static bool intersect_(
const Point& a,
const Point& b,
const Point& c,
1151 hasSegmentIs =
false;
1153 const auto ab = b - a;
1154 const auto ac = c - a;
1155 const auto qp = p - q;
1162 const auto denom = normal*qp;
1163 const auto abnorm2 = ab.two_norm2();
1164 const auto eps = eps_*abnorm2*qp.two_norm();
1166 if (abs(denom) < eps)
1168 const auto pa = a - p;
1169 const auto denom2 = normal*pa;
1170 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1175 using SegmentPolicy =
typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1176 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1177 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1178 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1179 using SegmentIntersectionType =
typename SegmentIntersectionAlgorithm::Intersection;
1180 SegmentIntersectionType segmentIs;
1182 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1183 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1184 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1186 hasSegmentIs =
true;
1192 for (
const auto& ip : {p, q})
1208 const auto ap = p - a;
1209 const auto t = (ap*normal)/denom;
1210 if (t < 0.0 - eps_)
return false;
1211 if (t > 1.0 + eps_)
return false;
1216 const auto v = (ac*e)/denom;
1217 if (v < -eps_ || v > 1.0 + eps_)
return false;
1218 const auto w = -(ab*e)/denom;
1219 if (w < -eps_ || v + w > 1.0 + eps_)
return false;
1223 const auto u = 1.0 - v - w;
1238template <
class Geometry1,
class Geometry2,
class Policy>
1251 template<
class P = Policy>
1262template <
class Geometry1,
class Geometry2,
class Policy>
1265 enum { dimworld = 3 };
1275 static constexpr ctype eps_ = 1.5e-7;
1276 using ReferenceElements =
typename Dune::ReferenceElements<ctype, dim1>;
1286 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1289 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1290 DUNE_THROW(Dune::NotImplemented,
"segment-segment intersection detection for point-like intersections");
1300 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1303 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1305 const auto& a = geo1.corner(0);
1306 const auto& b = geo1.corner(1);
1307 const auto& p = geo2.corner(0);
1308 const auto& q = geo2.corner(1);
1310 const auto ab = b-a;
1311 const auto pq = q-p;
1312 const auto abNorm = ab.two_norm();
1313 const auto pqNorm = pq.two_norm();
1316 const auto maxNorm = max(abNorm, pqNorm);
1317 const auto eps2 = eps_*maxNorm*maxNorm;
1321 if (!(abs(abs(ab*pq) - abNorm*pqNorm) < eps2))
1324 const auto ap = (p-a);
1325 const auto aq = (q-a);
1336 const auto abnorm2 = abNorm*abNorm;
1343 using std::min;
using std::max;
1344 tp = min(abnorm2, max(0.0, tp));
1345 tq = max(0.0, min(abnorm2, tq));
1347 if (abs(tp-tq) < eps2)
An axis-aligned bounding box volume hierarchy for dune grids.
Detect if a point intersects a geometry.
A function to compute the convex hull of a point cloud and a function to triangulate the polygon area...
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
bool intersectsPointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Check whether a point is intersectin a bounding box (dimworld == 3)
Definition: boundingboxtree.hh:304
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:114
Policy structure for point-like intersections.
Definition: geometryintersection.hh:44
static constexpr int dimIntersection
Definition: geometryintersection.hh:49
Point Intersection
Definition: geometryintersection.hh:50
ct ctype
Definition: geometryintersection.hh:45
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:46
static constexpr int dimWorld
Definition: geometryintersection.hh:48
Policy structure for segment-like intersections.
Definition: geometryintersection.hh:56
Segment Intersection
Definition: geometryintersection.hh:63
static constexpr int dimIntersection
Definition: geometryintersection.hh:62
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:58
static constexpr int dimWorld
Definition: geometryintersection.hh:61
std::array< Point, 2 > Segment
Definition: geometryintersection.hh:59
ct ctype
Definition: geometryintersection.hh:57
Policy structure for polygon-like intersections.
Definition: geometryintersection.hh:69
static constexpr int dimWorld
Definition: geometryintersection.hh:74
static constexpr int dimIntersection
Definition: geometryintersection.hh:75
ct ctype
Definition: geometryintersection.hh:70
Polygon Intersection
Definition: geometryintersection.hh:76
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:71
std::vector< Point > Polygon
Definition: geometryintersection.hh:72
Policy structure for polyhedron-like intersections.
Definition: geometryintersection.hh:82
std::vector< Point > Polyhedron
Definition: geometryintersection.hh:85
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:84
Polyhedron Intersection
Definition: geometryintersection.hh:89
static constexpr int dimIntersection
Definition: geometryintersection.hh:88
static constexpr int dimWorld
Definition: geometryintersection.hh:87
ct ctype
Definition: geometryintersection.hh:83
default policy chooser class
Definition: geometryintersection.hh:95
std::tuple_element_t< isDim, DefaultPolicies > type
Definition: geometryintersection.hh:109
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:217
typename Policy::ctype ctype
Definition: geometryintersection.hh:221
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:223
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:226
typename Policy::Point Point
Definition: geometryintersection.hh:222
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:251
typename Policy::ctype ctype
Definition: geometryintersection.hh:249
typename Policy::Point Point
Definition: geometryintersection.hh:250
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:260
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:380
typename Policy::ctype ctype
Definition: geometryintersection.hh:363
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:365
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:404
typename Policy::Point Point
Definition: geometryintersection.hh:364
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:491
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:511
typename Policy::ctype ctype
Definition: geometryintersection.hh:509
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:532
typename Policy::Point Point
Definition: geometryintersection.hh:510
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:657
typename Policy::Point Point
Definition: geometryintersection.hh:656
typename Policy::ctype ctype
Definition: geometryintersection.hh:655
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:676
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:786
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:806
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:829
typename Policy::Point Point
Definition: geometryintersection.hh:805
typename Policy::ctype ctype
Definition: geometryintersection.hh:804
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:978
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:998
typename Policy::Point Point
Definition: geometryintersection.hh:997
typename Policy::ctype ctype
Definition: geometryintersection.hh:996
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1015
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1041
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:1252
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1272
typename Policy::Point Point
Definition: geometryintersection.hh:1271
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:1287
typename Policy::ctype ctype
Definition: geometryintersection.hh:1270