3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
geometryintersection.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
23#ifndef DUMUX_GEOMETRY_INTERSECTION_HH
24#define DUMUX_GEOMETRY_INTERSECTION_HH
25
26#include <tuple>
27
28#include <dune/common/exceptions.hh>
29#include <dune/common/promotiontraits.hh>
30#include <dune/geometry/multilineargeometry.hh>
31
32#include <dumux/common/math.hh>
36
37namespace Dumux {
38namespace IntersectionPolicy {
39
41template<class ct, int dw>
43{
44 using ctype = ct;
45 using Point = Dune::FieldVector<ctype, dw>;
46
47 static constexpr int dimWorld = dw;
48 static constexpr int dimIntersection = 0;
50};
51
53template<class ct, int dw>
55{
56 using ctype = ct;
57 using Point = Dune::FieldVector<ctype, dw>;
58 using Segment = std::array<Point, 2>;
59
60 static constexpr int dimWorld = dw;
61 static constexpr int dimIntersection = 1;
63};
64
66template<class ct, int dw>
68{
69 using ctype = ct;
70 using Point = Dune::FieldVector<ctype, dw>;
71 using Polygon = std::vector<Point>;
72
73 static constexpr int dimWorld = dw;
74 static constexpr int dimIntersection = 2;
76};
77
79template<class ct, int dw>
81{
82 using ctype = ct;
83 using Point = Dune::FieldVector<ctype, dw>;
84 using Polyhedron = std::vector<Point>;
85
86 static constexpr int dimWorld = dw;
87 static constexpr int dimIntersection = 3;
89};
90
92template<class Geometry1, class Geometry2>
94{
95 static constexpr int dimworld = Geometry1::coorddimension;
96 static constexpr int isDim = std::min( int(Geometry1::mydimension), int(Geometry2::mydimension) );
97 static_assert(dimworld == int(Geometry2::coorddimension),
98 "Geometries must have the same coordinate dimension!");
99 static_assert(int(Geometry1::mydimension) <= 3 && int(Geometry2::mydimension) <= 3,
100 "Geometries must have dimension 3 or less.");
101 using ctype = typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
102
103 using DefaultPolicies = std::tuple<PointPolicy<ctype, dimworld>,
107public:
108 using type = std::tuple_element_t<isDim, DefaultPolicies>;
109};
110
112template<class Geometry1, class Geometry2>
114
115} // end namespace IntersectionPolicy
116
117namespace Detail {
118
132template< class Geo1, class Geo2, class ctype,
133 class GetFacetCornerIndices, class ComputeNormalFunction >
134bool computeSegmentIntersection(const Geo1& geo1, const Geo2& geo2, ctype baseEps,
135 ctype& tfirst, ctype& tlast,
136 const GetFacetCornerIndices& getFacetCornerIndices,
137 const ComputeNormalFunction& computeNormal)
138{
139 static_assert(int(Geo2::mydimension) == 1, "Geometry2 must be a segment");
140 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
141 "Dimension of Geometry1 must be higher than that of Geometry2");
142
143 const auto a = geo2.corner(0);
144 const auto b = geo2.corner(1);
145 const auto d = b - a;
146
147 // The initial interval is the whole segment.
148 // Afterwards we start clipping the interval
149 // at the intersections with the facets of geo1
150 tfirst = 0.0;
151 tlast = 1.0;
152
153 const auto facets = getFacetCornerIndices(geo1);
154 for (const auto& f : facets)
155 {
156 const auto n = computeNormal(f);
157
158 const auto c0 = geo1.corner(f[0]);
159 const ctype denom = n*d;
160 const ctype dist = n*(a-c0);
161
162 // use first edge of the facet to scale eps
163 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
164 const ctype eps = baseEps*edge1.two_norm();
165
166 // if denominator is zero the segment in parallel to the edge.
167 // If the distance is positive there is no intersection
168 using std::abs;
169 if (abs(denom) < eps)
170 {
171 if (dist > eps)
172 return false;
173 }
174 else // not parallel: compute line-line intersection
175 {
176 const ctype t = -dist / denom;
177 // if entering half space cut tfirst if t is larger
178 using std::signbit;
179 if (signbit(denom))
180 {
181 if (t > tfirst)
182 tfirst = t;
183 }
184 // if exiting half space cut tlast if t is smaller
185 else
186 {
187 if (t < tlast)
188 tlast = t;
189 }
190 // there is no intersection of the interval is empty
191 // use unscaled epsilon since t is in local coordinates
192 if (tfirst > tlast - baseEps)
193 return false;
194 }
195 }
196
197 // If we made it until here an intersection exists.
198 return true;
199}
200
201} // end namespace Detail
202
209template
210<class Geometry1, class Geometry2,
211 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
212 int dimworld = Geometry1::coorddimension,
213 int dim1 = Geometry1::mydimension,
214 int dim2 = Geometry2::mydimension>
216{
217 static constexpr int dimWorld = Policy::dimWorld;
218
219public:
220 using ctype = typename Policy::ctype;
221 using Point = typename Policy::Point;
222 using Intersection = typename Policy::Intersection;
223
225 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
226 {
227 static_assert(dimworld == Geometry2::coorddimension, "Can only intersect geometries of same coordinate dimension");
228 DUNE_THROW(Dune::NotImplemented, "Geometry intersection detection with intersection computation not implemented for dimworld = "
229 << dimworld << ", dim1 = " << dim1 << ", dim2 = " << dim2);
230 }
231};
232
237template <class Geometry1, class Geometry2, class Policy>
238class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 1>
239{
240 enum { dimworld = 2 };
241 enum { dim1 = 1 };
242 enum { dim2 = 1 };
243
244 // base epsilon for floating point comparisons
245 static constexpr typename Policy::ctype eps_ = 1.5e-7;
246
247public:
248 using ctype = typename Policy::ctype;
249 using Point = typename Policy::Point;
250 using Intersection = typename Policy::Intersection;
251
258 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
259 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
260 {
261 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
262
263 const auto v1 = geo1.corner(1) - geo1.corner(0);
264 const auto v2 = geo2.corner(1) - geo2.corner(0);
265 const auto ac = geo2.corner(0) - geo1.corner(0);
266
267 auto n2 = Point({-1.0*v2[1], v2[0]});
268 n2 /= n2.two_norm();
269
270 // compute distance of first corner on geo2 to line1
271 const auto dist12 = n2*ac;
272
273 // first check parallel segments
274 using std::abs;
275 using std::sqrt;
276
277 const auto v1Norm2 = v1.two_norm2();
278 const auto eps = eps_*sqrt(v1Norm2);
279 const auto eps2 = eps_*v1Norm2;
280
281 const auto sp = n2*v1;
282 if (abs(sp) < eps)
283 {
284 if (abs(dist12) > eps)
285 return false;
286
287 // point intersection can only be one of the corners
288 if ( (v1*v2) < 0.0 )
289 {
290 if ( ac.two_norm2() < eps2 )
291 { intersection = geo2.corner(0); return true; }
292
293 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
294 { intersection = geo2.corner(1); return true; }
295 }
296 else
297 {
298 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
299 { intersection = geo2.corner(1); return true; }
300
301 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
302 { intersection = geo2.corner(0); return true; }
303 }
304
305 // no intersection
306 return false;
307 }
308
309 // intersection point on v1 in local coords
310 const auto t1 = dist12 / sp;
311
312 // check if the local coords are valid
313 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
314 return false;
315
316 // compute global coordinates
317 auto isPoint = geo1.global(t1);
318
319 // check if point is in bounding box of v2
320 const auto c = geo2.corner(0);
321 const auto d = geo2.corner(1);
322
323 using std::min; using std::max;
324 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]) });
325
326 if ( intersectsPointBoundingBox(isPoint, bBox.data()) )
327 {
328 intersection = std::move(isPoint);
329 return true;
330 }
331
332 return false;
333 }
334
341 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
342 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
343 {
344 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
345
346 const auto& a = geo1.corner(0);
347 const auto& b = geo1.corner(1);
348 const auto ab = b - a;
349
350 const auto& p = geo2.corner(0);
351 const auto& q = geo2.corner(1);
352 const auto pq = q - p;
353 Dune::FieldVector<ctype, dimworld> n2{-pq[1], pq[0]};
354
355 using std::max;
356 const auto abNorm2 = ab.two_norm2();
357 const auto pqNorm2 = pq.two_norm2();
358 const auto eps2 = eps_*max(abNorm2, pqNorm2);
359
360 // non-parallel segments do not intersect in a segment.
361 using std::abs;
362 if (abs(n2*ab) > eps2)
363 return false;
364
365 // check if the segments lie on the same line
366 const auto ap = p - a;
367 if (abs(ap*n2) > eps2)
368 return false;
369
370 // compute scaled local coordinates of corner 1/2 of segment2 on segment1
371 auto t1 = ab*ap;
372 auto t2 = ab*(q - a);
373
374 using std::swap;
375 if (t1 > t2)
376 swap(t1, t2);
377
378 using std::clamp;
379 t1 = clamp(t1, 0.0, abNorm2);
380 t2 = clamp(t2, 0.0, abNorm2);
381
382 if (abs(t2-t1) < eps2)
383 return false;
384
385 intersection = Intersection({geo1.global(t1/abNorm2), geo1.global(t2/abNorm2)});
386 return true;
387 }
388};
389
394template <class Geometry1, class Geometry2, class Policy>
395class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 1>
396{
397 enum { dimworld = 2 };
398 enum { dim1 = 2 };
399 enum { dim2 = 1 };
400
401public:
402 using ctype = typename Policy::ctype;
403 using Point = typename Policy::Point;
404 using Intersection = typename Policy::Intersection;
405
406private:
407 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
408
409public:
417 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
418 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
419 {
420 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
421
422 ctype tfirst, tlast;
423 if (intersect_(geo1, geo2, tfirst, tlast))
424 {
425 // the intersection exists. Export the intersections geometry now:
426 // s(t) = a + t(b-a) in [tfirst, tlast]
427 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
428 return true;
429 }
430
431 return false;
432 }
433
441 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
442 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename P::Intersection& intersection)
443 {
444 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
445
446 ctype tfirst, tlast;
447 if (!intersect_(geo1, geo2, tfirst, tlast))
448 {
449 // if start & end point are same, the intersection is a point
450 using std::abs;
451 if ( tfirst > tlast - eps_ )
452 {
453 intersection = Intersection(geo2.global(tfirst));
454 return true;
455 }
456 }
457
458 return false;
459 }
460
461private:
465 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
466 {
467 // lambda to obtain the facet corners on geo1
468 auto getFacetCorners = [] (const Geometry1& geo1)
469 {
470 std::vector< std::array<int, 2> > facetCorners;
471 switch (geo1.corners())
472 {
473 case 4: // quadrilateral
474 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
475 break;
476 case 3: // triangle
477 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
478 break;
479 default:
480 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
481 << geo1.type() << " with "<< geo1.corners() << " corners.");
482 }
483
484 return facetCorners;
485 };
486
487 // lambda to obtain the normal on a facet
488 const auto center1 = geo1.center();
489 auto computeNormal = [&geo1, &center1] (const std::array<int, 2>& facetCorners)
490 {
491 const auto c0 = geo1.corner(facetCorners[0]);
492 const auto c1 = geo1.corner(facetCorners[1]);
493 const auto edge = c1 - c0;
494
495 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
496 n /= n.two_norm();
497
498 // orientation of the element is not known
499 // make sure normal points outwards of element
500 if ( n*(center1-c0) > 0.0 )
501 n *= -1.0;
502
503 return n;
504 };
505
506 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
507 }
508};
509
514template <class Geometry1, class Geometry2, class Policy>
515class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 2>
516: public GeometryIntersection<Geometry2, Geometry1, Policy, 2, 2, 1>
517{
519
520public:
528 template<class P = Policy>
529 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
530 {
531 return Base::intersection(geo2, geo1, intersection);
532 }
533};
534
539template <class Geometry1, class Geometry2, class Policy>
540class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 2>
541{
542 enum { dimworld = 2 };
543 enum { dim1 = 2 };
544 enum { dim2 = 2 };
545
546public:
547 using ctype = typename Policy::ctype;
548 using Point = typename Policy::Point;
549 using Intersection = typename Policy::Intersection;
550
551private:
552 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
553
554public:
567 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
568 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
569 {
570 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
571
572 // the candidate intersection points
573 std::vector<Point> points; points.reserve(6);
574
575 // add polygon1 corners that are inside polygon2
576 for (int i = 0; i < geo1.corners(); ++i)
577 if (intersectsPointGeometry(geo1.corner(i), geo2))
578 points.emplace_back(geo1.corner(i));
579
580 const auto numPoints1 = points.size();
581 if (numPoints1 != geo1.corners())
582 {
583 // add polygon2 corners that are inside polygon1
584 for (int i = 0; i < geo2.corners(); ++i)
585 if (intersectsPointGeometry(geo2.corner(i), geo1))
586 points.emplace_back(geo2.corner(i));
587
588 if (points.empty())
589 return false;
590
591 if (points.size() - numPoints1 != geo2.corners())
592 {
593 const auto refElement1 = referenceElement(geo1);
594 const auto refElement2 = referenceElement(geo2);
595
596 // add intersections of edges
597 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
599 for (int i = 0; i < refElement1.size(dim1-1); ++i)
600 {
601 const auto localEdgeGeom1 = refElement1.template geometry<dim1-1>(i);
602 const auto edge1 = SegGeometry( Dune::GeometryTypes::line,
603 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
604 geo1.global(localEdgeGeom1.corner(1))} ));
605
606 for (int j = 0; j < refElement2.size(dim2-1); ++j)
607 {
608 const auto localEdgeGeom2 = refElement2.template geometry<dim2-1>(j);
609 const auto edge2 = SegGeometry( Dune::GeometryTypes::line,
610 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
611 geo2.global(localEdgeGeom2.corner(1))} ));
612
614 typename EdgeTest::Intersection edgeIntersection;
615 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
616 points.emplace_back(edgeIntersection);
617 }
618 }
619 }
620 }
621
622 if (points.empty())
623 return false;
624
625 // remove duplicates
626 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
627 std::sort(points.begin(), points.end(), [&eps](const auto& a, const auto& b) -> bool
628 {
629 using std::abs;
630 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
631 });
632
633 auto removeIt = std::unique(points.begin(), points.end(), [&eps](const auto& a, const auto&b)
634 {
635 return (b-a).two_norm() < eps;
636 });
637
638 points.erase(removeIt, points.end());
639
640 // return false if we don't have at least three unique points
641 if (points.size() < 3)
642 return false;
643
644 // intersection polygon is convex hull of above points
645 intersection = grahamConvexHull<2>(points);
646 assert(!intersection.empty());
647 return true;
648 }
649
657 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
658 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
659 {
660 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
661 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for segment-like intersections");
662 }
663
671 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
672 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
673 {
674 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
675 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for touching points");
676 }
677};
678
683template <class Geometry1, class Geometry2, class Policy>
684class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 1>
685{
686 enum { dimworld = 3 };
687 enum { dim1 = 3 };
688 enum { dim2 = 1 };
689
690public:
691 using ctype = typename Policy::ctype;
692 using Point = typename Policy::Point;
693 using Intersection = typename Policy::Intersection;
694
695private:
696 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
697
698public:
710 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
711 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
712 {
713 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
714
715 ctype tfirst, tlast;
716 if (intersect_(geo1, geo2, tfirst, tlast))
717 {
718 // Intersection exists. Export the intersections geometry now:
719 // s(t) = a + t(b-a) in [tfirst, tlast]
720 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
721 return true;
722 }
723
724 return false;
725 }
726
738 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
739 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
740 {
741 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
742
743 ctype tfirst, tlast;
744 if (intersect_(geo1, geo2, tfirst, tlast))
745 {
746 // if start & end point are same, the intersection is a point
747 using std::abs;
748 if ( abs(tfirst - tlast) < eps_ )
749 {
750 intersection = Intersection(geo2.global(tfirst));
751 return true;
752 }
753 }
754
755 return false;
756 }
757
758private:
762 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
763 {
764 // lambda to obtain facet corners on geo1
765 auto getFacetCorners = [] (const Geometry1& geo1)
766 {
767 std::vector< std::vector<int> > facetCorners;
768 // sort facet corner so that normal n = (p1-p0)x(p2-p0) always points outwards
769 switch (geo1.corners())
770 {
771 // todo compute intersection geometries
772 case 8: // hexahedron
773 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
774 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
775 break;
776 case 6: // prism
777 facetCorners = {{0, 2, 1}, {3, 4, 5}, {0, 1, 3, 4}, {0, 3, 2, 5}, {1, 2, 4, 5}};
778 break;
779 case 5: // pyramid
780 facetCorners = {{0, 2, 1, 3}, {0, 1, 4}, {1, 3, 4}, {2, 4, 3}, {0, 4, 2}};
781 break;
782 case 4: // tetrahedron
783 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
784 break;
785 default:
786 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
787 << geo1.type() << ", "<< geo1.corners() << " corners.");
788 }
789
790 return facetCorners;
791 };
792
793 // lambda to obtain the normal on a facet
794 auto computeNormal = [&geo1] (const std::vector<int>& facetCorners)
795 {
796 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
797 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
798
799 auto n = crossProduct(v0, v1);
800 n /= n.two_norm();
801
802 return n;
803 };
804
805 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
806 }
807};
808
813template <class Geometry1, class Geometry2, class Policy>
814class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 3>
815: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 1>
816{
818public:
826 template<class P = Policy>
827 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
828 {
829 return Base::intersection(geo2, geo1, intersection);
830 }
831};
832
837template <class Geometry1, class Geometry2, class Policy>
838class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 2>
839{
840 enum { dimworld = 3 };
841 enum { dim1 = 2 };
842 enum { dim2 = 2 };
843
844public:
845 using ctype = typename Policy::ctype;
846 using Point = typename Policy::Point;
847 using Intersection = typename Policy::Intersection;
848
849private:
850 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
851 using ReferenceElementsGeo1 = typename Dune::ReferenceElements<ctype, dim1>;
852 using ReferenceElementsGeo2 = typename Dune::ReferenceElements<ctype, dim2>;
853
854public:
866 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
867 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
868 {
869 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
870
871 // if the geometries are not parallel, there cannot be a polygon intersection
872 const auto a1 = geo1.corner(1) - geo1.corner(0);
873 const auto b1 = geo1.corner(2) - geo1.corner(0);
874 const auto n1 = crossProduct(a1, b1);
875
876 const auto a2 = geo2.corner(1) - geo2.corner(0);
877 const auto b2 = geo2.corner(2) - geo2.corner(0);
878 const auto n2 = crossProduct(a2, b2);
879
880 // compute an epsilon scaled with larger edge length
881 const auto a1Norm2 = a1.two_norm2();
882 const auto b1Norm2 = b1.two_norm2();
883
884 using std::max;
885 const auto maxNorm2 = max(a1Norm2, b1Norm2);
886 const auto eps2 = maxNorm2*eps_;
887
888 // early exit if polygons are not parallel
889 using std::abs;
890 if (crossProduct(n1, n2).two_norm2() > eps2*maxNorm2)
891 return false;
892
893 // they are parallel, verify that they are on the same plane
894 auto d1 = geo2.corner(0) - geo1.corner(0);
895 if (d1.two_norm2() < eps2)
896 d1 = geo2.corner(1) - geo1.corner(0);
897
898 using std::sqrt;
899 const auto maxNorm = sqrt(maxNorm2);
900 const auto eps3 = maxNorm*eps2;
901
902 if (abs(d1*n2) > eps3)
903 return false;
904
905 // the candidate intersection points
906 std::vector<Point> points; points.reserve(8);
907
908 // add polygon1 corners that are inside polygon2
909 for (int i = 0; i < geo1.corners(); ++i)
910 if (intersectsPointGeometry(geo1.corner(i), geo2))
911 points.emplace_back(geo1.corner(i));
912
913 const auto numPoints1 = points.size();
914 const bool resultIsGeo1 = numPoints1 == geo1.corners();
915 if (!resultIsGeo1)
916 {
917 // add polygon2 corners that are inside polygon1
918 for (int i = 0; i < geo2.corners(); ++i)
919 if (intersectsPointGeometry(geo2.corner(i), geo1))
920 points.emplace_back(geo2.corner(i));
921
922 const bool resultIsGeo2 = (points.size() - numPoints1) == geo2.corners();
923 if (!resultIsGeo2)
924 {
925 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
926 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
927
928 // add intersections of edges
929 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
931 for (int i = 0; i < referenceElement1.size(dim1-1); ++i)
932 {
933 const auto localEdgeGeom1 = referenceElement1.template geometry<dim1-1>(i);
934 const auto edge1 = SegGeometry(
936 std::vector<Point>({
937 geo1.global(localEdgeGeom1.corner(0)),
938 geo1.global(localEdgeGeom1.corner(1))
939 })
940 );
941
942 for (int j = 0; j < referenceElement2.size(dim2-1); ++j)
943 {
944 const auto localEdgeGeom2 = referenceElement2.template geometry<dim2-1>(j);
945 const auto edge2 = SegGeometry(
947 std::vector<Point>({
948 geo2.global(localEdgeGeom2.corner(0)),
949 geo2.global(localEdgeGeom2.corner(1))
950 })
951 );
952
954 typename EdgeTest::Intersection edgeIntersection;
955 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
956 points.emplace_back(edgeIntersection);
957 }
958 }
959 }
960 }
961
962 if (points.empty())
963 return false;
964
965 // remove duplicates
966 const auto eps = maxNorm*eps_;
967 std::sort(points.begin(), points.end(), [eps] (const auto& a, const auto& b) -> bool
968 {
969 using std::abs;
970 return (abs(a[0]-b[0]) > eps ? a[0] < b[0]
971 : (abs(a[1]-b[1]) > eps ? a[1] < b[1]
972 : (a[2] < b[2])));
973 });
974
975 const auto squaredEps = eps*eps;
976 points.erase(std::unique(
977 points.begin(), points.end(),
978 [squaredEps] (const auto& a, const auto&b) { return (b-a).two_norm2() < squaredEps; }),
979 points.end()
980 );
981
982 // return false if we don't have at least three unique points
983 if (points.size() < 3)
984 return false;
985
986 // intersection polygon is convex hull of above points
987 intersection = grahamConvexHull<2>(points);
988 assert(!intersection.empty());
989 return true;
990 }
991
999 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1000 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1001 {
1002 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1003 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for segment-like intersections");
1004 }
1005
1013 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1014 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1015 {
1016 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1017 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for touching points");
1018 }
1019};
1020
1025template <class Geometry1, class Geometry2, class Policy>
1026class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 2>
1027{
1028 enum { dimworld = 3 };
1029 enum { dim1 = 3 };
1030 enum { dim2 = 2 };
1031
1032public:
1033 using ctype = typename Policy::ctype;
1034 using Point = typename Policy::Point;
1035 using Intersection = typename Policy::Intersection;
1036
1037private:
1038 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1039
1040public:
1055 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
1056 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1057 {
1058 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1059
1060 // the candidate intersection points
1061 std::vector<Point> points; points.reserve(10);
1062
1063 // add 3d geometry corners that are inside the 2d geometry
1064 for (int i = 0; i < geo1.corners(); ++i)
1065 if (intersectsPointGeometry(geo1.corner(i), geo2))
1066 points.emplace_back(geo1.corner(i));
1067
1068 // add 2d geometry corners that are inside the 3d geometry
1069 for (int i = 0; i < geo2.corners(); ++i)
1070 if (intersectsPointGeometry(geo2.corner(i), geo1))
1071 points.emplace_back(geo2.corner(i));
1072
1073 // get some geometry types
1074 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
1075 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
1076
1077 const auto refElement1 = referenceElement(geo1);
1078 const auto refElement2 = referenceElement(geo2);
1079
1080 // intersection policy for point-like intersections (used later)
1082
1083 // add intersection points of all polyhedron edges (codim dim-1) with the polygon
1084 for (int i = 0; i < refElement1.size(dim1-1); ++i)
1085 {
1086 const auto localEdgeGeom = refElement1.template geometry<dim1-1>(i);
1087 const auto p = geo1.global(localEdgeGeom.corner(0));
1088 const auto q = geo1.global(localEdgeGeom.corner(1));
1089 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
1090
1092 typename PolySegTest::Intersection polySegIntersection;
1093 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
1094 points.emplace_back(std::move(polySegIntersection));
1095 }
1096
1097 // add intersection points of all polygon faces (codim 1) with the polyhedron faces
1098 for (int i = 0; i < refElement1.size(1); ++i)
1099 {
1100 const auto faceGeo = [&]()
1101 {
1102 const auto localFaceGeo = refElement1.template geometry<1>(i);
1103 if (localFaceGeo.corners() == 4)
1104 {
1105 const auto a = geo1.global(localFaceGeo.corner(0));
1106 const auto b = geo1.global(localFaceGeo.corner(1));
1107 const auto c = geo1.global(localFaceGeo.corner(2));
1108 const auto d = geo1.global(localFaceGeo.corner(3));
1109
1110 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
1111 }
1112 else
1113 {
1114 const auto a = geo1.global(localFaceGeo.corner(0));
1115 const auto b = geo1.global(localFaceGeo.corner(1));
1116 const auto c = geo1.global(localFaceGeo.corner(2));
1117
1118 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
1119 }
1120 }();
1121
1122 for (int j = 0; j < refElement2.size(1); ++j)
1123 {
1124 const auto localEdgeGeom = refElement2.template geometry<1>(j);
1125 const auto p = geo2.global(localEdgeGeom.corner(0));
1126 const auto q = geo2.global(localEdgeGeom.corner(1));
1127
1128 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
1129
1131 typename PolySegTest::Intersection polySegIntersection;
1132 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
1133 points.emplace_back(std::move(polySegIntersection));
1134 }
1135 }
1136
1137 // return if no intersection points were found
1138 if (points.empty()) return false;
1139
1140 // remove duplicates
1141 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
1142 const auto notEqual = [eps] (auto a, auto b) { using std::abs; return abs(b-a) > eps; };
1143 std::sort(points.begin(), points.end(), [notEqual](const auto& a, const auto& b) -> bool
1144 {
1145 return (notEqual(a[0], b[0]) ? a[0] < b[0]
1146 : (notEqual(a[1], b[1]) ? a[1] < b[1]
1147 : (a[2] < b[2])));
1148 });
1149
1150 const auto squaredEps = eps*eps;
1151 points.erase(std::unique(
1152 points.begin(), points.end(),
1153 [squaredEps] (const auto& a, const auto&b) { return (b-a).two_norm2() < squaredEps; }),
1154 points.end()
1155 );
1156
1157 // return false if we don't have more than three unique points
1158 if (points.size() < 3) return false;
1159
1160 // intersection polygon is convex hull of above points
1161 intersection = grahamConvexHull<2>(points);
1162 assert(!intersection.empty());
1163
1164 return true;
1165 }
1166
1181 template<class P = Policy, std::enable_if_t<P::dimIntersection != 2, int> = 0>
1182 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1183 {
1184 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1185 DUNE_THROW(Dune::NotImplemented, "Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
1186 }
1187};
1188
1193template <class Geometry1, class Geometry2, class Policy>
1194class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 3>
1195: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 2>
1196{
1198public:
1206 template<class P = Policy>
1207 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
1208 {
1209 return Base::intersection(geo2, geo1, intersection);
1210 }
1211};
1212
1217template <class Geometry1, class Geometry2, class Policy>
1218class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 3>
1219{
1220 enum { dimworld = 3 };
1221 enum { dim1 = 3 };
1222 enum { dim2 = 3 };
1223
1224public:
1225 using ctype = typename Policy::ctype;
1226 using Point = typename Policy::Point;
1227 using Intersection = typename Policy::Intersection;
1228
1229private:
1230 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1231
1232public:
1244 template<class P = Policy, std::enable_if_t<P::dimIntersection == 3, int> = 0>
1245 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1246 {
1247 static_assert(
1248 int(dimworld) == int(Geometry2::coorddimension),
1249 "Can only collide geometries of same coordinate dimension"
1250 );
1251
1252 const auto refElement1 = referenceElement(geo1);
1253 const auto refElement2 = referenceElement(geo2);
1254
1255 // the candidate intersection points
1256 std::vector<Point> points;
1257 points.reserve(refElement1.size(2) + refElement2.size(2));
1258
1259 // add corners inside the other geometry
1260 const auto addPointIntersections = [&](const auto& g1, const auto& g2)
1261 {
1262 for (int i = 0; i < g1.corners(); ++i)
1263 if (const auto& c = g1.corner(i); intersectsPointGeometry(c, g2))
1264 points.emplace_back(c);
1265 };
1266
1267 addPointIntersections(geo1, geo2);
1268 addPointIntersections(geo2, geo1);
1269
1270 // get geometry types for the facets
1271 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
1272 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
1273
1274 // intersection policy for point-like intersections (used later)
1276
1277 // add intersection points of all edges with the faces of the other polyhedron
1278 const auto addEdgeIntersections = [&](const auto& g1, const auto& g2, const auto& ref1, const auto& ref2)
1279 {
1280 for (int i = 0; i < ref1.size(1); ++i)
1281 {
1282 const auto faceGeo = [&]()
1283 {
1284 const auto localFaceGeo = ref1.template geometry<1>(i);
1285 if (localFaceGeo.corners() == 4)
1286 {
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));
1290 const auto d = g1.global(localFaceGeo.corner(3));
1291
1292 return PolyhedronFaceGeometry(
1293 Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d}
1294 );
1295 }
1296 else
1297 {
1298 const auto a = g1.global(localFaceGeo.corner(0));
1299 const auto b = g1.global(localFaceGeo.corner(1));
1300 const auto c = g1.global(localFaceGeo.corner(2));
1301
1302 return PolyhedronFaceGeometry(
1303 Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c}
1304 );
1305 }
1306 }();
1307
1308 for (int j = 0; j < ref2.size(2); ++j)
1309 {
1310 const auto localEdgeGeom = ref2.template geometry<2>(j);
1311 const auto p = g2.global(localEdgeGeom.corner(0));
1312 const auto q = g2.global(localEdgeGeom.corner(1));
1313
1314 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
1315
1317 typename PolySegTest::Intersection polySegIntersection;
1318 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
1319 points.emplace_back(std::move(polySegIntersection));
1320 }
1321 }
1322 };
1323
1324 addEdgeIntersections(geo1, geo2, refElement1, refElement2);
1325 addEdgeIntersections(geo2, geo1, refElement2, refElement1);
1326
1327 // return if no intersection points were found
1328 if (points.empty())
1329 return false;
1330
1331 // remove duplicates
1332 const auto norm = (geo1.corner(0) - geo1.corner(1)).two_norm();
1333 const auto eps = norm*eps_;
1334 const auto notEqual = [eps] (auto a, auto b) { using std::abs; return abs(b-a) > eps; };
1335 std::sort(points.begin(), points.end(), [notEqual](const auto& a, const auto& b) -> bool
1336 {
1337 return (notEqual(a[0], b[0]) ? a[0] < b[0]
1338 : (notEqual(a[1], b[1]) ? a[1] < b[1]
1339 : (a[2] < b[2])));
1340 });
1341
1342 const auto squaredEps = eps*eps;
1343 points.erase(std::unique(
1344 points.begin(), points.end(),
1345 [squaredEps] (const auto& a, const auto&b) { return (b-a).two_norm2() < squaredEps; }),
1346 points.end()
1347 );
1348
1349 // return false if we don't have more than four unique points (dim+1)
1350 if (points.size() < 4)
1351 return false;
1352
1353 // check if the points are coplanar (then we don't have a 3-dimensional intersection)
1354 const bool coplanar = [&]
1355 {
1356 const auto p0 = points[0];
1357 const auto normal = crossProduct(points[1] - p0, points[2] - p0);
1358 // include the normal in eps (intead of norm*norm) since the normal can become very small
1359 // if the first three points are very close together
1360 const auto epsCoplanar = normal.two_norm()*norm*eps_;
1361 for (int i = 3; i < points.size(); ++i)
1362 {
1363 const auto ad = points[i] - p0;
1364 using std::abs;
1365 if (abs(normal*ad) > epsCoplanar)
1366 return false;
1367 }
1368
1369 return true;
1370 }();
1371
1372 if (coplanar)
1373 return false;
1374
1375 // we return the intersection as point cloud (that can be triangulated later)
1376 intersection = points;
1377
1378 return true;
1379 }
1380
1387 template<class P = Policy, std::enable_if_t<P::dimIntersection != 3, int> = 0>
1388 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1389 {
1390 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1391 DUNE_THROW(Dune::NotImplemented, "Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
1392 }
1393};
1394
1399template <class Geometry1, class Geometry2, class Policy>
1400class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 1>
1401{
1402 enum { dimworld = 3 };
1403 enum { dim1 = 2 };
1404 enum { dim2 = 1 };
1405
1406public:
1407 using ctype = typename Policy::ctype;
1408 using Point = typename Policy::Point;
1409 using Intersection = typename Policy::Intersection;
1410
1411private:
1412 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1413
1414public:
1424 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1425 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1426 {
1427 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1428
1429 ctype tfirst, tlast;
1430 if (intersect_(geo1, geo2, tfirst, tlast))
1431 {
1432 // the intersection exists. Export the intersections geometry now:
1433 // s(t) = a + t(b-a) in [tfirst, tlast]
1434 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
1435 return true;
1436 }
1437
1438 return false;
1439 }
1440
1450 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1451 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& is)
1452 {
1453 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1454
1455 const auto p = geo2.corner(0);
1456 const auto q = geo2.corner(1);
1457
1458 const auto a = geo1.corner(0);
1459 const auto b = geo1.corner(1);
1460 const auto c = geo1.corner(2);
1461
1462 if (geo1.corners() == 3)
1463 return intersect_<Policy>(a, b, c, p, q, is);
1464
1465 else if (geo1.corners() == 4)
1466 {
1467 Intersection is1, is2;
1468 bool hasSegment1, hasSegment2;
1469
1470 const auto d = geo1.corner(3);
1471 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1472 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1473
1474 if (hasSegment1 || hasSegment2)
1475 return false;
1476
1477 if (intersects1) { is = is1; return true; }
1478 if (intersects2) { is = is2; return true; }
1479
1480 return false;
1481 }
1482
1483 else
1484 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1485 << geo1.type() << ", "<< geo1.corners() << " corners.");
1486 }
1487
1488private:
1492 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1493 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
1494 {
1495 // lambda to obtain the facet corners on geo1
1496 auto getFacetCorners = [] (const Geometry1& geo1)
1497 {
1498 std::vector< std::array<int, 2> > facetCorners;
1499 switch (geo1.corners())
1500 {
1501 case 4: // quadrilateral
1502 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1503 break;
1504 case 3: // triangle
1505 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1506 break;
1507 default:
1508 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1509 << geo1.type() << " with "<< geo1.corners() << " corners.");
1510 }
1511
1512 return facetCorners;
1513 };
1514
1515 const auto center1 = geo1.center();
1516 const auto normal1 = crossProduct(geo1.corner(1) - geo1.corner(0),
1517 geo1.corner(2) - geo1.corner(0));
1518
1519 // lambda to obtain the normal on a facet
1520 auto computeNormal = [&center1, &normal1, &geo1] (const std::array<int, 2>& facetCorners)
1521 {
1522 const auto c0 = geo1.corner(facetCorners[0]);
1523 const auto c1 = geo1.corner(facetCorners[1]);
1524 const auto edge = c1 - c0;
1525
1526 Dune::FieldVector<ctype, dimworld> n = crossProduct(edge, normal1);
1527 n /= n.two_norm();
1528
1529 // orientation of the element is not known
1530 // make sure normal points outwards of element
1531 if ( n*(center1-c0) > 0.0 )
1532 n *= -1.0;
1533
1534 return n;
1535 };
1536
1537 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
1538 }
1539
1543 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1544 static bool intersect_(const Point& a, const Point& b, const Point& c,
1545 const Point& p, const Point& q,
1546 Intersection& is)
1547 {
1548 bool hasSegment;
1549 return intersect_(a, b, c, p, q, is, hasSegment);
1550 }
1551
1556 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1557 static bool intersect_(const Point& a, const Point& b, const Point& c,
1558 const Point& p, const Point& q,
1559 Intersection& is, bool& hasSegmentIs)
1560 {
1561 hasSegmentIs = false;
1562
1563 const auto ab = b - a;
1564 const auto ac = c - a;
1565 const auto qp = p - q;
1566
1567 // compute the triangle normal that defines the triangle plane
1568 const auto normal = crossProduct(ab, ac);
1569
1570 // compute the denominator
1571 // if denom is 0 the segment is parallel and we can return
1572 const auto denom = normal*qp;
1573 const auto abnorm2 = ab.two_norm2();
1574 const auto eps = eps_*abnorm2*qp.two_norm();
1575 using std::abs;
1576 if (abs(denom) < eps)
1577 {
1578 const auto pa = a - p;
1579 const auto denom2 = normal*pa;
1580 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1581 return false;
1582
1583 // if we get here, we are in-plane. Check if a
1584 // segment-like intersection with geometry 1 exists.
1585 using SegmentPolicy = typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1586 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1587 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1588 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1589 using SegmentIntersectionType = typename SegmentIntersectionAlgorithm::Intersection;
1590 SegmentIntersectionType segmentIs;
1591
1592 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1593 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1594 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1595 {
1596 hasSegmentIs = true;
1597 return false;
1598 }
1599
1600 // We are in-plane. A point-like
1601 // intersection can only be on the edges
1602 for (const auto& ip : {p, q})
1603 {
1604 if ( intersectsPointSimplex(ip, a, b)
1605 || intersectsPointSimplex(ip, b, c)
1606 || intersectsPointSimplex(ip, c, a) )
1607 {
1608 is = ip;
1609 return true;
1610 }
1611 }
1612
1613 return false;
1614 }
1615
1616 // compute intersection t value of pq with plane of triangle.
1617 // a segment intersects if and only if 0 <= t <= 1.
1618 const auto ap = p - a;
1619 const auto t = (ap*normal)/denom;
1620 if (t < 0.0 - eps_) return false;
1621 if (t > 1.0 + eps_) return false;
1622
1623 // compute the barycentric coordinates and check if the intersection point
1624 // is inside the bounds of the triangle
1625 const auto e = crossProduct(qp, ap);
1626 const auto v = (ac*e)/denom;
1627 if (v < -eps_ || v > 1.0 + eps_) return false;
1628 const auto w = -(ab*e)/denom;
1629 if (w < -eps_ || v + w > 1.0 + eps_) return false;
1630
1631 // Now we are sure there is an intersection points
1632 // Perform delayed division compute the last barycentric coordinate component
1633 const auto u = 1.0 - v - w;
1634
1635 Point ip(0.0);
1636 ip.axpy(u, a);
1637 ip.axpy(v, b);
1638 ip.axpy(w, c);
1639 is = ip;
1640 return true;
1641 }
1642};
1643
1648template <class Geometry1, class Geometry2, class Policy>
1649class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 2>
1650: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 2, 1>
1651{
1653public:
1661 template<class P = Policy>
1662 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
1663 {
1664 return Base::intersection(geo2, geo1, intersection);
1665 }
1666};
1667
1672template <class Geometry1, class Geometry2, class Policy>
1673class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 1>
1674{
1675 enum { dimworld = 3 };
1676 enum { dim1 = 1 };
1677 enum { dim2 = 1 };
1678
1679public:
1680 using ctype = typename Policy::ctype;
1681 using Point = typename Policy::Point;
1682 using Intersection = typename Policy::Intersection;
1683
1684private:
1685 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1686
1687public:
1694 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1695 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1696 {
1697 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1698
1699 const auto v1 = geo1.corner(1) - geo1.corner(0);
1700 const auto v2 = geo2.corner(1) - geo2.corner(0);
1701 const auto ac = geo2.corner(0) - geo1.corner(0);
1702
1703 const auto v1Norm2 = v1.two_norm2();
1704 const auto eps2 = eps_*v1Norm2;
1705
1706 const auto n = crossProduct(v1, v2);
1707
1708 // first check if segments are parallel
1709 using std::abs;
1710 if ( n.two_norm2() < eps2*v1Norm2 )
1711 {
1712 // check if they lie on the same line
1713 if (crossProduct(v1, ac).two_norm2() > eps2)
1714 return false;
1715
1716 // they lie on the same line,
1717 // if so, point intersection can only be one of the corners
1718 const auto sp = v1*v2;
1719 if ( sp < 0.0 )
1720 {
1721 if ( ac.two_norm2() < eps2 )
1722 { intersection = geo2.corner(0); return true; }
1723
1724 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
1725 { intersection = geo2.corner(1); return true; }
1726 }
1727 else
1728 {
1729 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
1730 { intersection = geo2.corner(1); return true; }
1731
1732 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
1733 { intersection = geo2.corner(0); return true; }
1734 }
1735
1736 // no intersection
1737 return false;
1738 }
1739
1740 // in-plane normal vector on v1
1741 const auto v1Normal = crossProduct(v1, n);
1742
1743 // intersection point on v2 in local coords
1744 const auto t2 = - 1.0*(ac*v1Normal) / (v2*v1Normal);
1745
1746 // check if the local coords are valid
1747 if (t2 < -1.0*eps_ || t2 > 1.0 + eps_)
1748 return false;
1749
1750 if (auto ip = geo2.global(t2); intersectsPointGeometry(ip, geo1))
1751 { intersection = std::move(ip); return true; }
1752
1753 return false;
1754 }
1755
1763 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1764 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1765 {
1766 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1767
1768 const auto& a = geo1.corner(0);
1769 const auto& b = geo1.corner(1);
1770 const auto ab = b-a;
1771
1772 const auto& p = geo2.corner(0);
1773 const auto& q = geo2.corner(1);
1774 const auto pq = q-p;
1775
1776 const auto abNorm2 = ab.two_norm2();
1777 const auto pqNorm2 = pq.two_norm2();
1778
1779 using std::max;
1780 const auto eps2 = eps_*max(abNorm2, pqNorm2);
1781
1782 // if the segment are not parallel there is no segment intersection
1783 using std::abs;
1784 if (crossProduct(ab, pq).two_norm2() > eps2*eps2)
1785 return false;
1786
1787 const auto ap = (p-a);
1788 const auto aq = (q-a);
1789
1790 // points have to be colinear
1791 if (crossProduct(ap, aq).two_norm2() > eps2*eps2)
1792 return false;
1793
1794 // scaled local coordinates
1795 // we save the division for the normalization for now
1796 // and do it in the very end if we find an intersection
1797 auto tp = ap*ab;
1798 auto tq = aq*ab;
1799
1800 // make sure they are sorted
1801 using std::swap;
1802 if (tp > tq)
1803 swap(tp, tq);
1804
1805 using std::clamp;
1806 tp = clamp(tp, 0.0, abNorm2);
1807 tq = clamp(tq, 0.0, abNorm2);
1808
1809 if (abs(tp-tq) < eps2)
1810 return false;
1811
1812 intersection = Intersection({geo1.global(tp/abNorm2), geo1.global(tq/abNorm2)});
1813 return true;
1814 }
1815};
1816
1817} // end namespace Dumux
1818
1819#endif
Define some often used mathematical functions.
Detect if a point intersects a geometry.
An axis-aligned bounding box volume hierarchy for dune grids.
A function to compute the convex hull of a point cloud.
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
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:36
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:36
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:134
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
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:304
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:113
constexpr Line line
Definition: couplingmanager1d3d_line.hh:43
Policy structure for point-like intersections.
Definition: geometryintersection.hh:43
static constexpr int dimIntersection
Definition: geometryintersection.hh:48
Point Intersection
Definition: geometryintersection.hh:49
ct ctype
Definition: geometryintersection.hh:44
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:45
static constexpr int dimWorld
Definition: geometryintersection.hh:47
Policy structure for segment-like intersections.
Definition: geometryintersection.hh:55
Segment Intersection
Definition: geometryintersection.hh:62
static constexpr int dimIntersection
Definition: geometryintersection.hh:61
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:57
static constexpr int dimWorld
Definition: geometryintersection.hh:60
std::array< Point, 2 > Segment
Definition: geometryintersection.hh:58
ct ctype
Definition: geometryintersection.hh:56
Policy structure for polygon-like intersections.
Definition: geometryintersection.hh:68
static constexpr int dimWorld
Definition: geometryintersection.hh:73
static constexpr int dimIntersection
Definition: geometryintersection.hh:74
ct ctype
Definition: geometryintersection.hh:69
Polygon Intersection
Definition: geometryintersection.hh:75
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:70
std::vector< Point > Polygon
Definition: geometryintersection.hh:71
Policy structure for polyhedron-like intersections.
Definition: geometryintersection.hh:81
std::vector< Point > Polyhedron
Definition: geometryintersection.hh:84
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:83
Polyhedron Intersection
Definition: geometryintersection.hh:88
static constexpr int dimIntersection
Definition: geometryintersection.hh:87
static constexpr int dimWorld
Definition: geometryintersection.hh:86
ct ctype
Definition: geometryintersection.hh:82
default policy chooser class
Definition: geometryintersection.hh:94
std::tuple_element_t< isDim, DefaultPolicies > type
Definition: geometryintersection.hh:108
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:216
typename Policy::ctype ctype
Definition: geometryintersection.hh:220
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:222
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:225
typename Policy::Point Point
Definition: geometryintersection.hh:221
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:250
typename Policy::ctype ctype
Definition: geometryintersection.hh:248
typename Policy::Point Point
Definition: geometryintersection.hh:249
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:259
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:418
typename Policy::ctype ctype
Definition: geometryintersection.hh:402
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:404
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:442
typename Policy::Point Point
Definition: geometryintersection.hh:403
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:529
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:549
typename Policy::ctype ctype
Definition: geometryintersection.hh:547
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:568
typename Policy::Point Point
Definition: geometryintersection.hh:548
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:693
typename Policy::Point Point
Definition: geometryintersection.hh:692
typename Policy::ctype ctype
Definition: geometryintersection.hh:691
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:711
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:827
typename Policy::Point Point
Definition: geometryintersection.hh:846
typename Policy::ctype ctype
Definition: geometryintersection.hh:845
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:847
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:867
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1035
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:1056
typename Policy::Point Point
Definition: geometryintersection.hh:1034
typename Policy::ctype ctype
Definition: geometryintersection.hh:1033
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:1207
typename Policy::Point Point
Definition: geometryintersection.hh:1226
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1227
typename Policy::ctype ctype
Definition: geometryintersection.hh:1225
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two convex polyhedra.
Definition: geometryintersection.hh:1245
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1409
typename Policy::Point Point
Definition: geometryintersection.hh:1408
typename Policy::ctype ctype
Definition: geometryintersection.hh:1407
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1425
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1451
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:1662
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1682
typename Policy::Point Point
Definition: geometryintersection.hh:1681
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:1695
typename Policy::ctype ctype
Definition: geometryintersection.hh:1680