25#ifndef DUMUX_GEOMETRY_INTERSECTION_HH
26#define DUMUX_GEOMETRY_INTERSECTION_HH
30#include <dune/common/exceptions.hh>
31#include <dune/common/promotiontraits.hh>
32#include <dune/geometry/multilineargeometry.hh>
33#include <dune/geometry/referenceelements.hh>
41namespace IntersectionPolicy {
44template<
class ct,
int dw>
48 using Point = Dune::FieldVector<ctype, dw>;
56template<
class ct,
int dw>
60 using Point = Dune::FieldVector<ctype, dw>;
69template<
class ct,
int dw>
73 using Point = Dune::FieldVector<ctype, dw>;
82template<
class ct,
int dw>
86 using Point = Dune::FieldVector<ctype, dw>;
95template<
class Geometry1,
class Geometry2>
98 static constexpr int dimworld = Geometry1::coorddimension;
99 static constexpr int isDim = std::min(
int(Geometry1::mydimension),
int(Geometry2::mydimension) );
100 static_assert(dimworld == int(Geometry2::coorddimension),
101 "Geometries must have the same coordinate dimension!");
102 static_assert(int(Geometry1::mydimension) <= 3 && int(Geometry2::mydimension) <= 3,
103 "Geometries must have dimension 3 or less.");
104 using ctype =
typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
106 using DefaultPolicies = std::tuple<PointPolicy<ctype, dimworld>,
111 using type = std::tuple_element_t<isDim, DefaultPolicies>;
115template<
class Geometry1,
class Geometry2>
135template<
class Geo1,
class Geo2,
class ctype,
136 class GetFacetCornerIndices,
class ComputeNormalFunction >
138 ctype& tfirst, ctype& tlast,
139 const GetFacetCornerIndices& getFacetCornerIndices,
140 const ComputeNormalFunction& computeNormal)
142 static_assert(int(Geo2::mydimension) == 1,
"Geometry2 must be a segment");
143 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
144 "Dimension of Geometry1 must be higher than that of Geometry2");
146 const auto a = geo2.corner(0);
147 const auto b = geo2.corner(1);
148 const auto d = b - a;
156 const auto facets = getFacetCornerIndices(geo1);
157 for (
const auto& f : facets)
159 const auto n = computeNormal(f);
161 const auto c0 = geo1.corner(f[0]);
162 const ctype denom = n*d;
163 const ctype dist = n*(a-c0);
166 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
167 const ctype eps = baseEps*edge1.two_norm();
172 if (abs(denom) < eps)
179 const ctype t = -dist / denom;
195 if (tfirst > tlast - baseEps)
213<
class Geometry1,
class Geometry2,
214 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
215 int dimworld = Geometry1::coorddimension,
216 int dim1 = Geometry1::mydimension,
217 int dim2 = Geometry2::mydimension>
220 static constexpr int dimWorld = Policy::dimWorld;
223 using ctype =
typename Policy::ctype;
224 using Point =
typename Policy::Point;
230 static_assert(dimworld == Geometry2::coorddimension,
"Can only intersect geometries of same coordinate dimension");
231 DUNE_THROW(Dune::NotImplemented,
"Geometry intersection detection with intersection computation not implemented for dimworld = "
232 << dimworld <<
", dim1 = " << dim1 <<
", dim2 = " << dim2);
240template <
class Geometry1,
class Geometry2,
class Policy>
243 enum { dimworld = 2 };
248 static constexpr typename Policy::ctype eps_ = 1.5e-7;
251 using ctype =
typename Policy::ctype;
252 using Point =
typename Policy::Point;
261 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
264 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
266 const auto v1 = geo1.corner(1) - geo1.corner(0);
267 const auto v2 = geo2.corner(1) - geo2.corner(0);
268 const auto ac = geo2.corner(0) - geo1.corner(0);
270 auto n2 =
Point({-1.0*v2[1], v2[0]});
274 const auto dist12 = n2*ac;
280 const auto v1Norm2 = v1.two_norm2();
281 const auto eps = eps_*sqrt(v1Norm2);
282 const auto eps2 = eps_*v1Norm2;
284 const auto sp = n2*v1;
287 if (abs(dist12) > eps)
293 if ( ac.two_norm2() < eps2 )
296 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
301 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
304 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
313 const auto t1 = dist12 / sp;
316 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
320 auto isPoint = geo1.global(t1);
323 const auto c = geo2.corner(0);
324 const auto d = geo2.corner(1);
326 using std::min;
using std::max;
327 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]) });
344 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
347 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
349 const auto& a = geo1.corner(0);
350 const auto& b = geo1.corner(1);
351 const auto ab = b - a;
353 const auto& p = geo2.corner(0);
354 const auto& q = geo2.corner(1);
355 const auto pq = q - p;
356 Dune::FieldVector<ctype, dimworld> n2{-pq[1], pq[0]};
359 const auto abNorm2 = ab.two_norm2();
360 const auto pqNorm2 = pq.two_norm2();
361 const auto eps2 = eps_*max(abNorm2, pqNorm2);
365 if (abs(n2*ab) > eps2)
369 const auto ap = p - a;
370 if (abs(ap*n2) > eps2)
375 auto t2 = ab*(q - a);
382 t1 = clamp(t1, 0.0, abNorm2);
383 t2 = clamp(t2, 0.0, abNorm2);
385 if (abs(t2-t1) < eps2)
397template <
class Geometry1,
class Geometry2,
class Policy>
400 enum { dimworld = 2 };
405 using ctype =
typename Policy::ctype;
406 using Point =
typename Policy::Point;
410 static constexpr ctype eps_ = 1.5e-7;
420 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
423 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
426 if (intersect_(geo1, geo2, tfirst, tlast))
444 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
447 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
450 if (!intersect_(geo1, geo2, tfirst, tlast))
454 if ( tfirst > tlast - eps_ )
468 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
471 auto getFacetCorners = [] (
const Geometry1& geo1)
473 std::vector< std::array<int, 2> > facetCorners;
474 switch (geo1.corners())
477 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
480 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
483 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
484 << geo1.type() <<
" with "<< geo1.corners() <<
" corners.");
491 const auto center1 = geo1.center();
492 auto computeNormal = [&geo1, ¢er1] (
const std::array<int, 2>& facetCorners)
494 const auto c0 = geo1.corner(facetCorners[0]);
495 const auto c1 = geo1.corner(facetCorners[1]);
496 const auto edge = c1 - c0;
498 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
503 if ( n*(center1-c0) > 0.0 )
517template <
class Geometry1,
class Geometry2,
class Policy>
531 template<
class P = Policy>
542template <
class Geometry1,
class Geometry2,
class Policy>
545 enum { dimworld = 2 };
550 using ctype =
typename Policy::ctype;
551 using Point =
typename Policy::Point;
555 static constexpr ctype eps_ = 1.5e-7;
570 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
573 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
576 std::vector<Point> points; points.reserve(6);
579 for (
int i = 0; i < geo1.corners(); ++i)
581 points.emplace_back(geo1.corner(i));
583 const auto numPoints1 = points.size();
584 if (numPoints1 != geo1.corners())
587 for (
int i = 0; i < geo2.corners(); ++i)
589 points.emplace_back(geo2.corner(i));
594 if (points.size() - numPoints1 != geo2.corners())
596 const auto refElement1 = referenceElement(geo1);
597 const auto refElement2 = referenceElement(geo2);
600 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
602 for (
int i = 0; i < refElement1.size(dim1-1); ++i)
604 const auto localEdgeGeom1 = refElement1.template geometry<dim1-1>(i);
606 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
607 geo1.global(localEdgeGeom1.corner(1))} ));
609 for (
int j = 0; j < refElement2.size(dim2-1); ++j)
611 const auto localEdgeGeom2 = refElement2.template geometry<dim2-1>(j);
613 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
614 geo2.global(localEdgeGeom2.corner(1))} ));
617 typename EdgeTest::Intersection edgeIntersection;
618 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
619 points.emplace_back(edgeIntersection);
629 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
630 std::sort(points.begin(), points.end(), [&eps](
const auto& a,
const auto& b) ->
bool
633 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
636 auto removeIt = std::unique(points.begin(), points.end(), [&eps](
const auto& a,
const auto&b)
638 return (b-a).two_norm() < eps;
641 points.erase(removeIt, points.end());
644 if (points.size() < 3)
660 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
663 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
664 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for segment-like intersections");
674 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
677 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
678 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for touching points");
686template <
class Geometry1,
class Geometry2,
class Policy>
689 enum { dimworld = 3 };
694 using ctype =
typename Policy::ctype;
695 using Point =
typename Policy::Point;
699 static constexpr ctype eps_ = 1.5e-7;
713 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
716 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
719 if (intersect_(geo1, geo2, tfirst, tlast))
741 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
744 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
747 if (intersect_(geo1, geo2, tfirst, tlast))
751 if ( abs(tfirst - tlast) < eps_ )
765 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
768 auto getFacetCorners = [] (
const Geometry1& geo1)
770 std::vector< std::vector<int> > facetCorners;
772 switch (geo1.corners())
776 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
777 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
780 facetCorners = {{0, 2, 1}, {3, 4, 5}, {0, 1, 3, 4}, {0, 3, 2, 5}, {1, 2, 4, 5}};
783 facetCorners = {{0, 2, 1, 3}, {0, 1, 4}, {1, 3, 4}, {2, 4, 3}, {0, 4, 2}};
786 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
789 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
790 << geo1.type() <<
", "<< geo1.corners() <<
" corners.");
797 auto computeNormal = [&geo1] (
const std::vector<int>& facetCorners)
799 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
800 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
816template <
class Geometry1,
class Geometry2,
class Policy>
829 template<
class P = Policy>
840template <
class Geometry1,
class Geometry2,
class Policy>
843 enum { dimworld = 3 };
848 using ctype =
typename Policy::ctype;
849 using Point =
typename Policy::Point;
853 static constexpr ctype eps_ = 1.5e-7;
867 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
870 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
873 const auto a1 = geo1.corner(1) - geo1.corner(0);
874 const auto b1 = geo1.corner(2) - geo1.corner(0);
877 const auto a2 = geo2.corner(1) - geo2.corner(0);
878 const auto b2 = geo2.corner(2) - geo2.corner(0);
882 const auto a1Norm2 = a1.two_norm2();
883 const auto b1Norm2 = b1.two_norm2();
886 const auto maxNorm2 = max(a1Norm2, b1Norm2);
887 const auto eps2 = maxNorm2*eps_;
895 auto d1 = geo2.corner(0) - geo1.corner(0);
896 if (d1.two_norm2() < eps2)
897 d1 = geo2.corner(1) - geo1.corner(0);
900 const auto maxNorm = sqrt(maxNorm2);
901 const auto eps3 = maxNorm*eps2;
903 if (abs(d1*n2) > eps3)
907 std::vector<Point> points; points.reserve(8);
910 for (
int i = 0; i < geo1.corners(); ++i)
912 points.emplace_back(geo1.corner(i));
914 const auto numPoints1 = points.size();
915 const bool resultIsGeo1 = numPoints1 == geo1.corners();
919 for (
int i = 0; i < geo2.corners(); ++i)
921 points.emplace_back(geo2.corner(i));
923 const bool resultIsGeo2 = (points.size() - numPoints1) == geo2.corners();
926 const auto referenceElement1 = referenceElement(geo1);
927 const auto referenceElement2 = referenceElement(geo2);
930 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
932 for (
int i = 0; i < referenceElement1.size(dim1-1); ++i)
934 const auto localEdgeGeom1 = referenceElement1.template geometry<dim1-1>(i);
935 const auto edge1 = SegGeometry(
938 geo1.global(localEdgeGeom1.corner(0)),
939 geo1.global(localEdgeGeom1.corner(1))
943 for (
int j = 0; j < referenceElement2.size(dim2-1); ++j)
945 const auto localEdgeGeom2 = referenceElement2.template geometry<dim2-1>(j);
946 const auto edge2 = SegGeometry(
949 geo2.global(localEdgeGeom2.corner(0)),
950 geo2.global(localEdgeGeom2.corner(1))
955 typename EdgeTest::Intersection edgeIntersection;
956 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
957 points.emplace_back(edgeIntersection);
967 const auto eps = maxNorm*eps_;
968 std::sort(points.begin(), points.end(), [eps] (
const auto& a,
const auto& b) ->
bool
971 return (abs(a[0]-b[0]) > eps ? a[0] < b[0]
972 : (abs(a[1]-b[1]) > eps ? a[1] < b[1]
976 const auto squaredEps = eps*eps;
977 points.erase(std::unique(
978 points.begin(), points.end(),
979 [squaredEps] (
const auto& a,
const auto&b) { return (b-a).two_norm2() < squaredEps; }),
984 if (points.size() < 3)
1000 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1003 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1004 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for segment-like intersections");
1014 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1017 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1018 DUNE_THROW(Dune::NotImplemented,
"Polygon-polygon intersection detection for touching points");
1026template <
class Geometry1,
class Geometry2,
class Policy>
1029 enum { dimworld = 3 };
1039 static constexpr ctype eps_ = 1.5e-7;
1056 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 2,
int> = 0>
1059 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1062 std::vector<Point> points; points.reserve(10);
1065 for (
int i = 0; i < geo1.corners(); ++i)
1067 points.emplace_back(geo1.corner(i));
1070 for (
int i = 0; i < geo2.corners(); ++i)
1072 points.emplace_back(geo2.corner(i));
1075 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
1076 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
1078 const auto refElement1 = referenceElement(geo1);
1079 const auto refElement2 = referenceElement(geo2);
1085 for (
int i = 0; i < refElement1.size(dim1-1); ++i)
1087 const auto localEdgeGeom = refElement1.template geometry<dim1-1>(i);
1088 const auto p = geo1.global(localEdgeGeom.corner(0));
1089 const auto q = geo1.global(localEdgeGeom.corner(1));
1093 typename PolySegTest::Intersection polySegIntersection;
1094 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
1095 points.emplace_back(std::move(polySegIntersection));
1099 for (
int i = 0; i < refElement1.size(1); ++i)
1101 const auto faceGeo = [&]()
1103 const auto localFaceGeo = refElement1.template geometry<1>(i);
1104 if (localFaceGeo.corners() == 4)
1106 const auto a = geo1.global(localFaceGeo.corner(0));
1107 const auto b = geo1.global(localFaceGeo.corner(1));
1108 const auto c = geo1.global(localFaceGeo.corner(2));
1109 const auto d = geo1.global(localFaceGeo.corner(3));
1111 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
1115 const auto a = geo1.global(localFaceGeo.corner(0));
1116 const auto b = geo1.global(localFaceGeo.corner(1));
1117 const auto c = geo1.global(localFaceGeo.corner(2));
1119 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
1123 for (
int j = 0; j < refElement2.size(1); ++j)
1125 const auto localEdgeGeom = refElement2.template geometry<1>(j);
1126 const auto p = geo2.global(localEdgeGeom.corner(0));
1127 const auto q = geo2.global(localEdgeGeom.corner(1));
1132 typename PolySegTest::Intersection polySegIntersection;
1133 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
1134 points.emplace_back(std::move(polySegIntersection));
1139 if (points.empty())
return false;
1142 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
1143 const auto notEqual = [eps] (
auto a,
auto b) {
using std::abs;
return abs(b-a) > eps; };
1144 std::sort(points.begin(), points.end(), [notEqual](
const auto& a,
const auto& b) ->
bool
1146 return (notEqual(a[0], b[0]) ? a[0] < b[0]
1147 : (notEqual(a[1], b[1]) ? a[1] < b[1]
1151 const auto squaredEps = eps*eps;
1152 points.erase(std::unique(
1153 points.begin(), points.end(),
1154 [squaredEps] (
const auto& a,
const auto&b) { return (b-a).two_norm2() < squaredEps; }),
1159 if (points.size() < 3)
return false;
1182 template<
class P = Policy, std::enable_if_t<P::dimIntersection != 2,
int> = 0>
1185 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1186 DUNE_THROW(Dune::NotImplemented,
"Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
1194template <
class Geometry1,
class Geometry2,
class Policy>
1207 template<
class P = Policy>
1218template <
class Geometry1,
class Geometry2,
class Policy>
1221 enum { dimworld = 3 };
1231 static constexpr ctype eps_ = 1.5e-7;
1245 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 3,
int> = 0>
1249 int(dimworld) == int(Geometry2::coorddimension),
1250 "Can only collide geometries of same coordinate dimension"
1253 const auto refElement1 = referenceElement(geo1);
1254 const auto refElement2 = referenceElement(geo2);
1257 std::vector<Point> points;
1258 points.reserve(refElement1.size(2) + refElement2.size(2));
1261 const auto addPointIntersections = [&](
const auto& g1,
const auto& g2)
1263 for (
int i = 0; i < g1.corners(); ++i)
1265 points.emplace_back(c);
1268 addPointIntersections(geo1, geo2);
1269 addPointIntersections(geo2, geo1);
1272 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
1273 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
1279 const auto addEdgeIntersections = [&](
const auto& g1,
const auto& g2,
const auto& ref1,
const auto& ref2)
1281 for (
int i = 0; i < ref1.size(1); ++i)
1283 const auto faceGeo = [&]()
1285 const auto localFaceGeo = ref1.template geometry<1>(i);
1286 if (localFaceGeo.corners() == 4)
1288 const auto a = g1.global(localFaceGeo.corner(0));
1289 const auto b = g1.global(localFaceGeo.corner(1));
1290 const auto c = g1.global(localFaceGeo.corner(2));
1291 const auto d = g1.global(localFaceGeo.corner(3));
1293 return PolyhedronFaceGeometry(
1294 Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d}
1299 const auto a = g1.global(localFaceGeo.corner(0));
1300 const auto b = g1.global(localFaceGeo.corner(1));
1301 const auto c = g1.global(localFaceGeo.corner(2));
1303 return PolyhedronFaceGeometry(
1304 Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c}
1309 for (
int j = 0; j < ref2.size(2); ++j)
1311 const auto localEdgeGeom = ref2.template geometry<2>(j);
1312 const auto p = g2.global(localEdgeGeom.corner(0));
1313 const auto q = g2.global(localEdgeGeom.corner(1));
1318 typename PolySegTest::Intersection polySegIntersection;
1319 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
1320 points.emplace_back(std::move(polySegIntersection));
1325 addEdgeIntersections(geo1, geo2, refElement1, refElement2);
1326 addEdgeIntersections(geo2, geo1, refElement2, refElement1);
1333 const auto norm = (geo1.corner(0) - geo1.corner(1)).two_norm();
1334 const auto eps = norm*eps_;
1335 const auto notEqual = [eps] (
auto a,
auto b) {
using std::abs;
return abs(b-a) > eps; };
1336 std::sort(points.begin(), points.end(), [notEqual](
const auto& a,
const auto& b) ->
bool
1338 return (notEqual(a[0], b[0]) ? a[0] < b[0]
1339 : (notEqual(a[1], b[1]) ? a[1] < b[1]
1343 const auto squaredEps = eps*eps;
1344 points.erase(std::unique(
1345 points.begin(), points.end(),
1346 [squaredEps] (
const auto& a,
const auto&b) { return (b-a).two_norm2() < squaredEps; }),
1351 if (points.size() < 4)
1355 const bool coplanar = [&]
1357 const auto p0 = points[0];
1361 const auto epsCoplanar =
normal.two_norm()*norm*eps_;
1362 for (
int i = 3; i < points.size(); ++i)
1364 const auto ad = points[i] - p0;
1366 if (abs(
normal*ad) > epsCoplanar)
1388 template<
class P = Policy, std::enable_if_t<P::dimIntersection != 3,
int> = 0>
1391 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1392 DUNE_THROW(Dune::NotImplemented,
"Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
1400template <
class Geometry1,
class Geometry2,
class Policy>
1403 enum { dimworld = 3 };
1413 static constexpr ctype eps_ = 1.5e-7;
1425 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1428 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1430 ctype tfirst, tlast;
1431 if (intersect_(geo1, geo2, tfirst, tlast))
1451 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1454 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1456 const auto p = geo2.corner(0);
1457 const auto q = geo2.corner(1);
1459 const auto a = geo1.corner(0);
1460 const auto b = geo1.corner(1);
1461 const auto c = geo1.corner(2);
1463 if (geo1.corners() == 3)
1464 return intersect_<Policy>(a, b, c, p, q, is);
1466 else if (geo1.corners() == 4)
1469 bool hasSegment1, hasSegment2;
1471 const auto d = geo1.corner(3);
1472 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1473 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1475 if (hasSegment1 || hasSegment2)
1478 if (intersects1) { is = is1;
return true; }
1479 if (intersects2) { is = is2;
return true; }
1485 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
1486 << geo1.type() <<
", "<< geo1.corners() <<
" corners.");
1493 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1494 static bool intersect_(
const Geometry1& geo1,
const Geometry2& geo2,
ctype& tfirst,
ctype& tlast)
1497 auto getFacetCorners = [] (
const Geometry1& geo1)
1499 std::vector< std::array<int, 2> > facetCorners;
1500 switch (geo1.corners())
1503 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1506 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1509 DUNE_THROW(Dune::NotImplemented,
"Collision of segment and geometry of type "
1510 << geo1.type() <<
" with "<< geo1.corners() <<
" corners.");
1513 return facetCorners;
1516 const auto center1 = geo1.center();
1517 const auto normal1 =
crossProduct(geo1.corner(1) - geo1.corner(0),
1518 geo1.corner(2) - geo1.corner(0));
1521 auto computeNormal = [¢er1, &normal1, &geo1] (
const std::array<int, 2>& facetCorners)
1523 const auto c0 = geo1.corner(facetCorners[0]);
1524 const auto c1 = geo1.corner(facetCorners[1]);
1525 const auto edge = c1 - c0;
1527 Dune::FieldVector<ctype, dimworld> n =
crossProduct(edge, normal1);
1532 if ( n*(center1-c0) > 0.0 )
1544 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1545 static bool intersect_(
const Point& a,
const Point& b,
const Point& c,
1550 return intersect_(a, b, c, p, q, is, hasSegment);
1557 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1558 static bool intersect_(
const Point& a,
const Point& b,
const Point& c,
1562 hasSegmentIs =
false;
1564 const auto ab = b - a;
1565 const auto ac = c - a;
1566 const auto qp = p - q;
1573 const auto denom =
normal*qp;
1574 const auto abnorm2 = ab.two_norm2();
1575 const auto eps = eps_*abnorm2*qp.two_norm();
1577 if (abs(denom) < eps)
1579 const auto pa = a - p;
1580 const auto denom2 =
normal*pa;
1581 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1586 using SegmentPolicy =
typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1587 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1588 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1589 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1590 using SegmentIntersectionType =
typename SegmentIntersectionAlgorithm::Intersection;
1591 SegmentIntersectionType segmentIs;
1593 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1594 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1595 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1597 hasSegmentIs =
true;
1603 for (
const auto& ip : {p, q})
1619 const auto ap = p - a;
1620 const auto t = (ap*
normal)/denom;
1621 if (t < 0.0 - eps_)
return false;
1622 if (t > 1.0 + eps_)
return false;
1627 const auto v = (ac*e)/denom;
1628 if (v < -eps_ || v > 1.0 + eps_)
return false;
1629 const auto w = -(ab*e)/denom;
1630 if (w < -eps_ || v + w > 1.0 + eps_)
return false;
1634 const auto u = 1.0 - v - w;
1649template <
class Geometry1,
class Geometry2,
class Policy>
1662 template<
class P = Policy>
1673template <
class Geometry1,
class Geometry2,
class Policy>
1676 enum { dimworld = 3 };
1686 static constexpr ctype eps_ = 1.5e-7;
1695 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 0,
int> = 0>
1698 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1700 const auto v1 = geo1.corner(1) - geo1.corner(0);
1701 const auto v2 = geo2.corner(1) - geo2.corner(0);
1702 const auto ac = geo2.corner(0) - geo1.corner(0);
1704 const auto v1Norm2 = v1.two_norm2();
1705 const auto eps2 = eps_*v1Norm2;
1711 if ( n.two_norm2() < eps2*v1Norm2 )
1719 const auto sp = v1*v2;
1722 if ( ac.two_norm2() < eps2 )
1725 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
1730 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
1733 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
1745 const auto t2 = - 1.0*(ac*v1Normal) / (v2*v1Normal);
1748 if (t2 < -1.0*eps_ || t2 > 1.0 + eps_)
1764 template<
class P = Policy, std::enable_if_t<P::dimIntersection == 1,
int> = 0>
1767 static_assert(int(dimworld) == int(Geometry2::coorddimension),
"Can only collide geometries of same coordinate dimension");
1769 const auto& a = geo1.corner(0);
1770 const auto& b = geo1.corner(1);
1771 const auto ab = b-a;
1773 const auto& p = geo2.corner(0);
1774 const auto& q = geo2.corner(1);
1775 const auto pq = q-p;
1777 const auto abNorm2 = ab.two_norm2();
1778 const auto pqNorm2 = pq.two_norm2();
1781 const auto eps2 = eps_*max(abNorm2, pqNorm2);
1788 const auto ap = (p-a);
1789 const auto aq = (q-a);
1807 tp = clamp(tp, 0.0, abNorm2);
1808 tq = clamp(tq, 0.0, abNorm2);
1810 if (abs(tp-tq) < eps2)
Define some often used mathematical functions.
An axis-aligned bounding box volume hierarchy for dune grids.
A function to compute the convex hull of a point cloud.
Detect if a point intersects a geometry.
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:40
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.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 the tetrahedron (p0, p1, p2, p3) (dimworld is 3)
Definition: intersectspointsimplex.hh:38
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:137
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:654
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
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: boundingboxtree.hh:306
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:116
constexpr Line line
Definition: couplingmanager1d3d_line.hh:43
Policy structure for point-like intersections.
Definition: geometryintersection.hh:46
static constexpr int dimIntersection
Definition: geometryintersection.hh:51
Point Intersection
Definition: geometryintersection.hh:52
ct ctype
Definition: geometryintersection.hh:47
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:48
static constexpr int dimWorld
Definition: geometryintersection.hh:50
Policy structure for segment-like intersections.
Definition: geometryintersection.hh:58
Segment Intersection
Definition: geometryintersection.hh:65
static constexpr int dimIntersection
Definition: geometryintersection.hh:64
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:60
static constexpr int dimWorld
Definition: geometryintersection.hh:63
std::array< Point, 2 > Segment
Definition: geometryintersection.hh:61
ct ctype
Definition: geometryintersection.hh:59
Policy structure for polygon-like intersections.
Definition: geometryintersection.hh:71
static constexpr int dimWorld
Definition: geometryintersection.hh:76
static constexpr int dimIntersection
Definition: geometryintersection.hh:77
ct ctype
Definition: geometryintersection.hh:72
Polygon Intersection
Definition: geometryintersection.hh:78
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:73
std::vector< Point > Polygon
Definition: geometryintersection.hh:74
Policy structure for polyhedron-like intersections.
Definition: geometryintersection.hh:84
std::vector< Point > Polyhedron
Definition: geometryintersection.hh:87
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:86
Polyhedron Intersection
Definition: geometryintersection.hh:91
static constexpr int dimIntersection
Definition: geometryintersection.hh:90
static constexpr int dimWorld
Definition: geometryintersection.hh:89
ct ctype
Definition: geometryintersection.hh:85
default policy chooser class
Definition: geometryintersection.hh:97
std::tuple_element_t< isDim, DefaultPolicies > type
Definition: geometryintersection.hh:111
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:219
typename Policy::ctype ctype
Definition: geometryintersection.hh:223
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:225
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:228
typename Policy::Point Point
Definition: geometryintersection.hh:224
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:253
typename Policy::ctype ctype
Definition: geometryintersection.hh:251
typename Policy::Point Point
Definition: geometryintersection.hh:252
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:262
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:421
typename Policy::ctype ctype
Definition: geometryintersection.hh:405
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:407
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:445
typename Policy::Point Point
Definition: geometryintersection.hh:406
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:532
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:552
typename Policy::ctype ctype
Definition: geometryintersection.hh:550
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:571
typename Policy::Point Point
Definition: geometryintersection.hh:551
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:696
typename Policy::Point Point
Definition: geometryintersection.hh:695
typename Policy::ctype ctype
Definition: geometryintersection.hh:694
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:714
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:830
typename Policy::Point Point
Definition: geometryintersection.hh:849
typename Policy::ctype ctype
Definition: geometryintersection.hh:848
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:850
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:868
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1036
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:1057
typename Policy::Point Point
Definition: geometryintersection.hh:1035
typename Policy::ctype ctype
Definition: geometryintersection.hh:1034
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:1208
typename Policy::Point Point
Definition: geometryintersection.hh:1227
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1228
typename Policy::ctype ctype
Definition: geometryintersection.hh:1226
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two convex polyhedra.
Definition: geometryintersection.hh:1246
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1410
typename Policy::Point Point
Definition: geometryintersection.hh:1409
typename Policy::ctype ctype
Definition: geometryintersection.hh:1408
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1426
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1452
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:1663
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1683
typename Policy::Point Point
Definition: geometryintersection.hh:1682
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:1696
typename Policy::ctype ctype
Definition: geometryintersection.hh:1681