13#ifndef DUMUX_GEOMETRY_INTERSECTION_HH
14#define DUMUX_GEOMETRY_INTERSECTION_HH
18#include <dune/common/exceptions.hh>
19#include <dune/common/promotiontraits.hh>
20#include <dune/geometry/multilineargeometry.hh>
21#include <dune/geometry/referenceelements.hh>
29namespace IntersectionPolicy {
32template<
class ct,
int dw>
36 using Point = Dune::FieldVector<ctype, dw>;
44template<
class ct,
int dw>
48 using Point = Dune::FieldVector<ctype, dw>;
57template<
class ct,
int dw>
61 using Point = Dune::FieldVector<ctype, dw>;
70template<
class ct,
int dw>
74 using Point = Dune::FieldVector<ctype, dw>;
83template<
class Geometry1,
class Geometry2>
86 static constexpr int dimworld = Geometry1::coorddimension;
87 static constexpr int isDim = std::min(
int(Geometry1::mydimension),
int(Geometry2::mydimension) );
88 static_assert(dimworld == int(Geometry2::coorddimension),
89 "Geometries must have the same coordinate dimension!");
90 static_assert(int(Geometry1::mydimension) <= 3 && int(Geometry2::mydimension) <= 3,
91 "Geometries must have dimension 3 or less.");
92 using ctype =
typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
94 using DefaultPolicies = std::tuple<PointPolicy<ctype, dimworld>,
99 using type = std::tuple_element_t<isDim, DefaultPolicies>;
103template<
class Geometry1,
class Geometry2>
123template<
class Geo1,
class Geo2,
class ctype,
124 class GetFacetCornerIndices,
class ComputeNormalFunction >
126 ctype& tfirst, ctype& tlast,
127 const GetFacetCornerIndices& getFacetCornerIndices,
128 const ComputeNormalFunction& computeNormal)
130 static_assert(int(Geo2::mydimension) == 1,
"Geometry2 must be a segment");
131 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
132 "Dimension of Geometry1 must be higher than that of Geometry2");
134 const auto a = geo2.corner(0);
135 const auto b = geo2.corner(1);
136 const auto d = b - a;
144 const auto facets = getFacetCornerIndices(geo1);
145 for (
const auto& f : facets)
147 const auto n = computeNormal(f);
149 const auto c0 = geo1.corner(f[0]);
150 const ctype denom = n*d;
151 const ctype dist = n*(a-c0);
154 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
155 const ctype eps = baseEps*edge1.two_norm();
160 if (abs(denom) < eps)
167 const ctype t = -dist / denom;
183 if (tfirst > tlast - baseEps)
201<
class Geometry1,
class Geometry2,
202 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
203 int dimworld = Geometry1::coorddimension,
204 int dim1 = Geometry1::mydimension,
205 int dim2 = Geometry2::mydimension>
208 static constexpr int dimWorld = Policy::dimWorld;
211 using ctype =
typename Policy::ctype;
212 using Point =
typename Policy::Point;
218 static_assert(dimworld == Geometry2::coorddimension,
"Can only intersect geometries of same coordinate dimension");
219 DUNE_THROW(Dune::NotImplemented,
"Geometry intersection detection with intersection computation not implemented for dimworld = "
220 << dimworld <<
", dim1 = " << dim1 <<
", dim2 = " << dim2);
228template <
class Geometry1,
class Geometry2,
class Policy>
231 enum { dimworld = 2 };
236 static constexpr typename Policy::ctype eps_ = 1.5e-7;
239 using ctype =
typename Policy::ctype;
240 using Point =
typename Policy::Point;
249 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
252 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
254 const auto v1 = geo1.corner(1) - geo1.corner(0);
255 const auto v2 = geo2.corner(1) - geo2.corner(0);
256 const auto ac = geo2.corner(0) - geo1.corner(0);
258 auto n2 =
Point({-1.0*v2[1], v2[0]});
262 const auto dist12 = n2*ac;
268 const auto v1Norm2 = v1.two_norm2();
269 const auto eps = eps_*sqrt(v1Norm2);
270 const auto eps2 = eps_*v1Norm2;
272 const auto sp = n2*v1;
275 if (abs(dist12) > eps)
281 if ( ac.two_norm2() < eps2 )
284 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
289 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
292 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
301 const auto t1 = dist12 / sp;
304 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
308 auto isPoint = geo1.global(t1);
311 const auto c = geo2.corner(0);
312 const auto d = geo2.corner(1);
314 using std::min;
using std::max;
315 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]) });
332 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
335 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
337 const auto& a = geo1.corner(0);
338 const auto& b = geo1.corner(1);
339 const auto ab = b - a;
341 const auto& p = geo2.corner(0);
342 const auto& q = geo2.corner(1);
343 const auto pq = q - p;
344 Dune::FieldVector<ctype, dimworld> n2{-pq[1], pq[0]};
347 const auto abNorm2 = ab.two_norm2();
348 const auto pqNorm2 = pq.two_norm2();
349 const auto eps2 = eps_*max(abNorm2, pqNorm2);
353 if (abs(n2*ab) > eps2)
357 const auto ap = p - a;
358 if (abs(ap*n2) > eps2)
363 auto t2 = ab*(q - a);
370 t1 = clamp(t1, 0.0, abNorm2);
371 t2 = clamp(t2, 0.0, abNorm2);
373 if (abs(t2-t1) < eps2)
385template <
class Geometry1,
class Geometry2,
class Policy>
388 enum { dimworld = 2 };
393 using ctype =
typename Policy::ctype;
394 using Point =
typename Policy::Point;
398 static constexpr ctype eps_ = 1.5e-7;
408 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
411 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
414 if (intersect_(geo1, geo2, tfirst, tlast))
432 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
435 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
438 if (!intersect_(geo1, geo2, tfirst, tlast))
442 if ( tfirst > tlast - eps_ )
456 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
459 auto getFacetCorners = [] (
const Geometry1& geo1)
461 std::vector< std::array<int, 2> > facetCorners;
462 switch (geo1.corners())
465 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
468 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
471 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
472 << geo1.type() <<
" with "<< geo1.corners() <<
" corners.");
479 const auto center1 = geo1.center();
480 auto computeNormal = [&geo1, ¢er1] (
const std::array<int, 2>& facetCorners)
482 const auto c0 = geo1.corner(facetCorners[0]);
483 const auto c1 = geo1.corner(facetCorners[1]);
484 const auto edge = c1 - c0;
486 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
491 if ( n*(center1-c0) > 0.0 )
505template <
class Geometry1,
class Geometry2,
class Policy>
519 template<
class P = Policy>
530template <
class Geometry1,
class Geometry2,
class Policy>
533 enum { dimworld = 2 };
538 using ctype =
typename Policy::ctype;
539 using Point =
typename Policy::Point;
543 static constexpr ctype eps_ = 1.5e-7;
558 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
561 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
564 std::vector<Point> points; points.reserve(6);
567 for (
int i = 0; i < geo1.corners(); ++i)
569 points.emplace_back(geo1.corner(i));
571 const auto numPoints1 = points.size();
572 if (numPoints1 != geo1.corners())
575 for (
int i = 0; i < geo2.corners(); ++i)
577 points.emplace_back(geo2.corner(i));
582 if (points.size() - numPoints1 != geo2.corners())
584 const auto refElement1 = referenceElement(geo1);
585 const auto refElement2 = referenceElement(geo2);
588 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
590 for (
int i = 0; i < refElement1.size(dim1-1); ++i)
592 const auto localEdgeGeom1 = refElement1.template geometry<dim1-1>(i);
594 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
595 geo1.global(localEdgeGeom1.corner(1))} ));
597 for (
int j = 0; j < refElement2.size(dim2-1); ++j)
599 const auto localEdgeGeom2 = refElement2.template geometry<dim2-1>(j);
601 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
602 geo2.global(localEdgeGeom2.corner(1))} ));
605 typename EdgeTest::Intersection edgeIntersection;
606 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
607 points.emplace_back(edgeIntersection);
617 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
618 std::sort(points.begin(), points.end(), [&eps](
const auto& a,
const auto& b) ->
bool
621 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
624 auto removeIt = std::unique(points.begin(), points.end(), [&eps](
const auto& a,
const auto&b)
626 return (b-a).two_norm() < eps;
629 points.erase(removeIt, points.end());
632 if (points.size() < 3)
648 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
651 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
652 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for segment-like intersections");
662 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
665 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
666 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for touching points");
674template <
class Geometry1,
class Geometry2,
class Policy>
677 enum { dimworld = 3 };
682 using ctype =
typename Policy::ctype;
683 using Point =
typename Policy::Point;
687 static constexpr ctype eps_ = 1.5e-7;
701 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
704 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
707 if (intersect_(geo1, geo2, tfirst, tlast))
729 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
732 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
735 if (intersect_(geo1, geo2, tfirst, tlast))
739 if ( abs(tfirst - tlast) < eps_ )
753 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
756 auto getFacetCorners = [] (
const Geometry1& geo1)
758 std::vector< std::vector<int> > facetCorners;
760 switch (geo1.corners())
764 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
765 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
768 facetCorners = {{0, 2, 1}, {3, 4, 5}, {0, 1, 3, 4}, {0, 3, 2, 5}, {1, 2, 4, 5}};
771 facetCorners = {{0, 2, 1, 3}, {0, 1, 4}, {1, 3, 4}, {2, 4, 3}, {0, 4, 2}};
774 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
777 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
778 << geo1.type() <<
", "<< geo1.corners() <<
" corners.");
785 auto computeNormal = [&geo1] (
const std::vector<int>& facetCorners)
787 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
788 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
804template <
class Geometry1,
class Geometry2,
class Policy>
817 template<
class P = Policy>
828template <
class Geometry1,
class Geometry2,
class Policy>
831 enum { dimworld = 3 };
836 using ctype =
typename Policy::ctype;
837 using Point =
typename Policy::Point;
841 static constexpr ctype eps_ = 1.5e-7;
855 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
858 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
861 const auto a1 = geo1.corner(1) - geo1.corner(0);
862 const auto b1 = geo1.corner(2) - geo1.corner(0);
865 const auto a2 = geo2.corner(1) - geo2.corner(0);
866 const auto b2 = geo2.corner(2) - geo2.corner(0);
870 const auto a1Norm2 = a1.two_norm2();
871 const auto b1Norm2 = b1.two_norm2();
874 const auto maxNorm2 = max(a1Norm2, b1Norm2);
875 const auto eps2 = maxNorm2*eps_;
883 auto d1 = geo2.corner(0) - geo1.corner(0);
884 if (d1.two_norm2() < eps2)
885 d1 = geo2.corner(1) - geo1.corner(0);
888 const auto maxNorm = sqrt(maxNorm2);
889 const auto eps3 = maxNorm*eps2;
891 if (abs(d1*n2) > eps3)
895 std::vector<Point> points; points.reserve(8);
898 for (
int i = 0; i < geo1.corners(); ++i)
900 points.emplace_back(geo1.corner(i));
902 const auto numPoints1 = points.size();
903 const bool resultIsGeo1 = numPoints1 == geo1.corners();
907 for (
int i = 0; i < geo2.corners(); ++i)
909 points.emplace_back(geo2.corner(i));
911 const bool resultIsGeo2 = (points.size() - numPoints1) == geo2.corners();
914 const auto referenceElement1 = referenceElement(geo1);
915 const auto referenceElement2 = referenceElement(geo2);
918 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
920 for (
int i = 0; i < referenceElement1.size(dim1-1); ++i)
922 const auto localEdgeGeom1 = referenceElement1.template geometry<dim1-1>(i);
923 const auto edge1 = SegGeometry(
926 geo1.global(localEdgeGeom1.corner(0)),
927 geo1.global(localEdgeGeom1.corner(1))
931 for (
int j = 0; j < referenceElement2.size(dim2-1); ++j)
933 const auto localEdgeGeom2 = referenceElement2.template geometry<dim2-1>(j);
934 const auto edge2 = SegGeometry(
937 geo2.global(localEdgeGeom2.corner(0)),
938 geo2.global(localEdgeGeom2.corner(1))
943 typename EdgeTest::Intersection edgeIntersection;
944 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
945 points.emplace_back(edgeIntersection);
955 const auto eps = maxNorm*eps_;
956 std::sort(points.begin(), points.end(), [eps] (
const auto& a,
const auto& b) ->
bool
959 return (abs(a[0]-b[0]) > eps ? a[0] < b[0]
960 : (abs(a[1]-b[1]) > eps ? a[1] < b[1]
964 const auto squaredEps = eps*eps;
965 points.erase(std::unique(
966 points.begin(), points.end(),
967 [squaredEps] (
const auto& a,
const auto&b) { return (b-a).two_norm2() < squaredEps; }),
972 if (points.size() < 3)
988 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
991 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
992 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for segment-like intersections");
1002 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1005 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1006 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for touching points");
1014template <
class Geometry1,
class Geometry2,
class Policy>
1017 enum { dimworld = 3 };
1027 static constexpr ctype eps_ = 1.5e-7;
1044 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
1047 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1050 std::vector<Point> points; points.reserve(10);
1053 for (
int i = 0; i < geo1.corners(); ++i)
1055 points.emplace_back(geo1.corner(i));
1058 for (
int i = 0; i < geo2.corners(); ++i)
1060 points.emplace_back(geo2.corner(i));
1063 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
1064 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
1066 const auto refElement1 = referenceElement(geo1);
1067 const auto refElement2 = referenceElement(geo2);
1073 for (
int i = 0; i < refElement1.size(dim1-1); ++i)
1075 const auto localEdgeGeom = refElement1.template geometry<dim1-1>(i);
1076 const auto p = geo1.global(localEdgeGeom.corner(0));
1077 const auto q = geo1.global(localEdgeGeom.corner(1));
1081 typename PolySegTest::Intersection polySegIntersection;
1082 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
1083 points.emplace_back(std::move(polySegIntersection));
1087 for (
int i = 0; i < refElement1.size(1); ++i)
1089 const auto faceGeo = [&]()
1091 const auto localFaceGeo = refElement1.template geometry<1>(i);
1092 if (localFaceGeo.corners() == 4)
1094 const auto a = geo1.global(localFaceGeo.corner(0));
1095 const auto b = geo1.global(localFaceGeo.corner(1));
1096 const auto c = geo1.global(localFaceGeo.corner(2));
1097 const auto d = geo1.global(localFaceGeo.corner(3));
1099 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
1103 const auto a = geo1.global(localFaceGeo.corner(0));
1104 const auto b = geo1.global(localFaceGeo.corner(1));
1105 const auto c = geo1.global(localFaceGeo.corner(2));
1107 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
1111 for (
int j = 0; j < refElement2.size(1); ++j)
1113 const auto localEdgeGeom = refElement2.template geometry<1>(j);
1114 const auto p = geo2.global(localEdgeGeom.corner(0));
1115 const auto q = geo2.global(localEdgeGeom.corner(1));
1120 typename PolySegTest::Intersection polySegIntersection;
1121 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
1122 points.emplace_back(std::move(polySegIntersection));
1127 if (points.empty())
return false;
1130 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
1131 const auto notEqual = [eps] (
auto a,
auto b) {
using std::abs;
return abs(b-a) > eps; };
1132 std::sort(points.begin(), points.end(), [notEqual](
const auto& a,
const auto& b) ->
bool
1134 return (notEqual(a[0], b[0]) ? a[0] < b[0]
1135 : (notEqual(a[1], b[1]) ? a[1] < b[1]
1139 const auto squaredEps = eps*eps;
1140 points.erase(std::unique(
1141 points.begin(), points.end(),
1142 [squaredEps] (
const auto& a,
const auto&b) { return (b-a).two_norm2() < squaredEps; }),
1147 if (points.size() < 3)
return false;
1170 template<
class P = Policy, std::enable_if_t<P::dimIntersection != 2,
int> = 0>
1173 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1174 DUNE_THROW(Dune::NotImplemented,
"Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
1182template <
class Geometry1,
class Geometry2,
class Policy>
1195 template<
class P = Policy>
1206template <
class Geometry1,
class Geometry2,
class Policy>
1209 enum { dimworld = 3 };
1219 static constexpr ctype eps_ = 1.5e-7;
1233 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 3,
int> = 0>
1237 int(dimworld) == int(Geometry2::coorddimension),
1238 "Can only collide geometries of same coordinate dimension"
1241 const auto refElement1 = referenceElement(geo1);
1242 const auto refElement2 = referenceElement(geo2);
1245 std::vector<Point> points;
1246 points.reserve(refElement1.size(2) + refElement2.size(2));
1249 const auto addPointIntersections = [&](
const auto& g1,
const auto& g2)
1251 for (
int i = 0; i < g1.corners(); ++i)
1253 points.emplace_back(c);
1256 addPointIntersections(geo1, geo2);
1257 addPointIntersections(geo2, geo1);
1260 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
1261 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
1267 const auto addEdgeIntersections = [&](
const auto& g1,
const auto& g2,
const auto& ref1,
const auto& ref2)
1269 for (
int i = 0; i < ref1.size(1); ++i)
1271 const auto faceGeo = [&]()
1273 const auto localFaceGeo = ref1.template geometry<1>(i);
1274 if (localFaceGeo.corners() == 4)
1276 const auto a = g1.global(localFaceGeo.corner(0));
1277 const auto b = g1.global(localFaceGeo.corner(1));
1278 const auto c = g1.global(localFaceGeo.corner(2));
1279 const auto d = g1.global(localFaceGeo.corner(3));
1281 return PolyhedronFaceGeometry(
1282 Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d}
1287 const auto a = g1.global(localFaceGeo.corner(0));
1288 const auto b = g1.global(localFaceGeo.corner(1));
1289 const auto c = g1.global(localFaceGeo.corner(2));
1291 return PolyhedronFaceGeometry(
1292 Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c}
1297 for (
int j = 0; j < ref2.size(2); ++j)
1299 const auto localEdgeGeom = ref2.template geometry<2>(j);
1300 const auto p = g2.global(localEdgeGeom.corner(0));
1301 const auto q = g2.global(localEdgeGeom.corner(1));
1306 typename PolySegTest::Intersection polySegIntersection;
1307 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
1308 points.emplace_back(std::move(polySegIntersection));
1313 addEdgeIntersections(geo1, geo2, refElement1, refElement2);
1314 addEdgeIntersections(geo2, geo1, refElement2, refElement1);
1321 const auto norm = (geo1.corner(0) - geo1.corner(1)).two_norm();
1322 const auto eps = norm*eps_;
1323 const auto notEqual = [eps] (
auto a,
auto b) {
using std::abs;
return abs(b-a) > eps; };
1324 std::sort(points.begin(), points.end(), [notEqual](
const auto& a,
const auto& b) ->
bool
1326 return (notEqual(a[0], b[0]) ? a[0] < b[0]
1327 : (notEqual(a[1], b[1]) ? a[1] < b[1]
1331 const auto squaredEps = eps*eps;
1332 points.erase(std::unique(
1333 points.begin(), points.end(),
1334 [squaredEps] (
const auto& a,
const auto&b) { return (b-a).two_norm2() < squaredEps; }),
1339 if (points.size() < 4)
1343 const bool coplanar = [&]
1345 const auto p0 = points[0];
1349 const auto epsCoplanar =
normal.two_norm()*norm*eps_;
1350 for (
int i = 3; i < points.size(); ++i)
1352 const auto ad = points[i] - p0;
1354 if (abs(
normal*ad) > epsCoplanar)
1376 template<
class P = Policy, std::enable_if_t<P::dimIntersection != 3,
int> = 0>
1379 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1380 DUNE_THROW(Dune::NotImplemented,
"Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
1388template <
class Geometry1,
class Geometry2,
class Policy>
1391 enum { dimworld = 3 };
1401 static constexpr ctype eps_ = 1.5e-7;
1413 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1416 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1418 ctype tfirst, tlast;
1419 if (intersect_(geo1, geo2, tfirst, tlast))
1439 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1442 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1444 const auto p = geo2.corner(0);
1445 const auto q = geo2.corner(1);
1447 const auto a = geo1.corner(0);
1448 const auto b = geo1.corner(1);
1449 const auto c = geo1.corner(2);
1451 if (geo1.corners() == 3)
1452 return intersect_<Policy>(a, b, c, p, q, is);
1454 else if (geo1.corners() == 4)
1457 bool hasSegment1, hasSegment2;
1459 const auto d = geo1.corner(3);
1460 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1461 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1463 if (hasSegment1 || hasSegment2)
1466 if (intersects1) { is = is1;
return true; }
1467 if (intersects2) { is = is2;
return true; }
1473 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
1474 << geo1.type() <<
", "<< geo1.corners() <<
" corners.");
1481 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1482 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
1485 auto getFacetCorners = [] (
const Geometry1& geo1)
1487 std::vector< std::array<int, 2> > facetCorners;
1488 switch (geo1.corners())
1491 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1494 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1497 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
1498 << geo1.type() <<
" with "<< geo1.corners() <<
" corners.");
1501 return facetCorners;
1504 const auto center1 = geo1.center();
1505 const auto normal1 =
crossProduct(geo1.corner(1) - geo1.corner(0),
1506 geo1.corner(2) - geo1.corner(0));
1509 auto computeNormal = [¢er1, &normal1, &geo1] (
const std::array<int, 2>& facetCorners)
1511 const auto c0 = geo1.corner(facetCorners[0]);
1512 const auto c1 = geo1.corner(facetCorners[1]);
1513 const auto edge = c1 - c0;
1515 Dune::FieldVector<ctype, dimworld> n =
crossProduct(edge, normal1);
1520 if ( n*(center1-c0) > 0.0 )
1532 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1533 static bool intersect_(
const Point& a,
const Point& b,
const Point& c,
1538 return intersect_(a, b, c, p, q, is, hasSegment);
1545 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1546 static bool intersect_(
const Point& a,
const Point& b,
const Point& c,
1550 hasSegmentIs =
false;
1552 const auto ab = b - a;
1553 const auto ac = c - a;
1554 const auto qp = p - q;
1561 const auto denom =
normal*qp;
1562 const auto abnorm2 = ab.two_norm2();
1563 const auto eps = eps_*abnorm2*qp.two_norm();
1565 if (abs(denom) < eps)
1567 const auto pa = a - p;
1568 const auto denom2 =
normal*pa;
1569 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1574 using SegmentPolicy =
typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1575 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1576 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1577 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1578 using SegmentIntersectionType =
typename SegmentIntersectionAlgorithm::Intersection;
1579 SegmentIntersectionType segmentIs;
1581 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1582 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1583 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1585 hasSegmentIs =
true;
1591 for (
const auto& ip : {p, q})
1607 const auto ap = p - a;
1608 const auto t = (ap*
normal)/denom;
1609 if (t < 0.0 - eps_)
return false;
1610 if (t > 1.0 + eps_)
return false;
1615 const auto v = (ac*e)/denom;
1616 if (v < -eps_ || v > 1.0 + eps_)
return false;
1617 const auto w = -(ab*e)/denom;
1618 if (w < -eps_ || v + w > 1.0 + eps_)
return false;
1622 const auto u = 1.0 - v - w;
1637template <
class Geometry1,
class Geometry2,
class Policy>
1650 template<
class P = Policy>
1661template <
class Geometry1,
class Geometry2,
class Policy>
1664 enum { dimworld = 3 };
1674 static constexpr ctype eps_ = 1.5e-7;
1683 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1686 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1688 const auto v1 = geo1.corner(1) - geo1.corner(0);
1689 const auto v2 = geo2.corner(1) - geo2.corner(0);
1690 const auto ac = geo2.corner(0) - geo1.corner(0);
1692 const auto v1Norm2 = v1.two_norm2();
1693 const auto eps2 = eps_*v1Norm2;
1699 if ( n.two_norm2() < eps2*v1Norm2 )
1707 const auto sp = v1*v2;
1710 if ( ac.two_norm2() < eps2 )
1713 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
1718 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
1721 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
1733 const auto t2 = - 1.0*(ac*v1Normal) / (v2*v1Normal);
1736 if (t2 < -1.0*eps_ || t2 > 1.0 + eps_)
1752 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1755 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1757 const auto& a = geo1.corner(0);
1758 const auto& b = geo1.corner(1);
1759 const auto ab = b-a;
1761 const auto& p = geo2.corner(0);
1762 const auto& q = geo2.corner(1);
1763 const auto pq = q-p;
1765 const auto abNorm2 = ab.two_norm2();
1766 const auto pqNorm2 = pq.two_norm2();
1769 const auto eps2 = eps_*max(abNorm2, pqNorm2);
1776 const auto ap = (p-a);
1777 const auto aq = (q-a);
1795 tp = clamp(tp, 0.0, abNorm2);
1796 tq = clamp(tq, 0.0, abNorm2);
1798 if (abs(tp-tq) < eps2)
An axis-aligned bounding box volume hierarchy for dune grids.
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:241
typename Policy::ctype ctype
Definition: geometryintersection.hh:239
typename Policy::Point Point
Definition: geometryintersection.hh:240
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:250
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:520
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:409
typename Policy::ctype ctype
Definition: geometryintersection.hh:393
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:395
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:433
typename Policy::Point Point
Definition: geometryintersection.hh:394
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:540
typename Policy::ctype ctype
Definition: geometryintersection.hh:538
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:559
typename Policy::Point Point
Definition: geometryintersection.hh:539
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1671
typename Policy::Point Point
Definition: geometryintersection.hh:1670
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:1684
typename Policy::ctype ctype
Definition: geometryintersection.hh:1669
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:1651
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:818
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1398
typename Policy::Point Point
Definition: geometryintersection.hh:1397
typename Policy::ctype ctype
Definition: geometryintersection.hh:1396
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1414
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1440
typename Policy::Point Point
Definition: geometryintersection.hh:837
typename Policy::ctype ctype
Definition: geometryintersection.hh:836
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:838
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:856
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:1196
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:684
typename Policy::Point Point
Definition: geometryintersection.hh:683
typename Policy::ctype ctype
Definition: geometryintersection.hh:682
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:702
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1024
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:1045
typename Policy::Point Point
Definition: geometryintersection.hh:1023
typename Policy::ctype ctype
Definition: geometryintersection.hh:1022
typename Policy::Point Point
Definition: geometryintersection.hh:1215
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1216
typename Policy::ctype ctype
Definition: geometryintersection.hh:1214
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two convex polyhedra.
Definition: geometryintersection.hh:1234
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:207
typename Policy::ctype ctype
Definition: geometryintersection.hh:211
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:213
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:216
typename Policy::Point Point
Definition: geometryintersection.hh:212
default policy chooser class
Definition: geometryintersection.hh:85
std::tuple_element_t< isDim, DefaultPolicies > type
Definition: geometryintersection.hh:99
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
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:28
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
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 the tetrahedron (p0, p1, p2, p3) (dimworld is 3)
Definition: intersectspointsimplex.hh:26
bool computeSegmentIntersection(const Geo1 &geo1, const Geo2 &geo2, ctype baseEps, ctype &tfirst, ctype &tlast, const GetFacetCornerIndices &getFacetCornerIndices, const ComputeNormalFunction &computeNormal)
Algorithm to find segment-like intersections of a polygon/polyhedron with a segment....
Definition: geometryintersection.hh:125
Detect if a point intersects a geometry.
Define some often used mathematical functions.
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:104
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:294
Policy structure for point-like intersections.
Definition: geometryintersection.hh:34
static constexpr int dimIntersection
Definition: geometryintersection.hh:39
Point Intersection
Definition: geometryintersection.hh:40
ct ctype
Definition: geometryintersection.hh:35
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:36
static constexpr int dimWorld
Definition: geometryintersection.hh:38
Policy structure for polygon-like intersections.
Definition: geometryintersection.hh:59
static constexpr int dimWorld
Definition: geometryintersection.hh:64
static constexpr int dimIntersection
Definition: geometryintersection.hh:65
ct ctype
Definition: geometryintersection.hh:60
Polygon Intersection
Definition: geometryintersection.hh:66
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:61
std::vector< Point > Polygon
Definition: geometryintersection.hh:62
Policy structure for polyhedron-like intersections.
Definition: geometryintersection.hh:72
std::vector< Point > Polyhedron
Definition: geometryintersection.hh:75
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:74
Polyhedron Intersection
Definition: geometryintersection.hh:79
static constexpr int dimIntersection
Definition: geometryintersection.hh:78
static constexpr int dimWorld
Definition: geometryintersection.hh:77
ct ctype
Definition: geometryintersection.hh:73
Policy structure for segment-like intersections.
Definition: geometryintersection.hh:46
Segment Intersection
Definition: geometryintersection.hh:53
static constexpr int dimIntersection
Definition: geometryintersection.hh:52
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:48
static constexpr int dimWorld
Definition: geometryintersection.hh:51
std::array< Point, 2 > Segment
Definition: geometryintersection.hh:49
ct ctype
Definition: geometryintersection.hh:47