version 3.10-dev
geometryintersection.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_GEOMETRY_INTERSECTION_HH
14#define DUMUX_GEOMETRY_INTERSECTION_HH
15
16#include <tuple>
17
18#include <dune/common/exceptions.hh>
19#include <dune/common/promotiontraits.hh>
20#include <dune/geometry/multilineargeometry.hh>
21#include <dune/geometry/referenceelements.hh>
22
23#include <dumux/common/math.hh>
27
28namespace Dumux {
29namespace IntersectionPolicy {
30
32template<class ct, int dw>
34{
35 using ctype = ct;
36 using Point = Dune::FieldVector<ctype, dw>;
37
38 static constexpr int dimWorld = dw;
39 static constexpr int dimIntersection = 0;
41};
42
44template<class ct, int dw>
46{
47 using ctype = ct;
48 using Point = Dune::FieldVector<ctype, dw>;
49 using Segment = std::array<Point, 2>;
50
51 static constexpr int dimWorld = dw;
52 static constexpr int dimIntersection = 1;
54};
55
57template<class ct, int dw>
59{
60 using ctype = ct;
61 using Point = Dune::FieldVector<ctype, dw>;
62 using Polygon = std::vector<Point>;
63
64 static constexpr int dimWorld = dw;
65 static constexpr int dimIntersection = 2;
67};
68
70template<class ct, int dw>
72{
73 using ctype = ct;
74 using Point = Dune::FieldVector<ctype, dw>;
75 using Polyhedron = std::vector<Point>;
76
77 static constexpr int dimWorld = dw;
78 static constexpr int dimIntersection = 3;
80};
81
83template<class Geometry1, class Geometry2>
85{
86 static constexpr int dimworld = Geometry1::coorddimension;
87 static constexpr int isDim = std::min( int(Geometry1::mydimension), int(Geometry2::mydimension) );
88 static_assert(dimworld == int(Geometry2::coorddimension),
89 "Geometries must have the same coordinate dimension!");
90 static_assert(int(Geometry1::mydimension) <= 3 && int(Geometry2::mydimension) <= 3,
91 "Geometries must have dimension 3 or less.");
92 using ctype = typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
93
94 using DefaultPolicies = std::tuple<PointPolicy<ctype, dimworld>,
98public:
99 using type = std::tuple_element_t<isDim, DefaultPolicies>;
100};
101
103template<class Geometry1, class Geometry2>
105
106} // end namespace IntersectionPolicy
107
108namespace Detail {
109
123template< class Geo1, class Geo2, class ctype,
124 class GetFacetCornerIndices, class ComputeNormalFunction >
125bool computeSegmentIntersection(const Geo1& geo1, const Geo2& geo2, ctype baseEps,
126 ctype& tfirst, ctype& tlast,
127 const GetFacetCornerIndices& getFacetCornerIndices,
128 const ComputeNormalFunction& computeNormal)
129{
130 static_assert(int(Geo2::mydimension) == 1, "Geometry2 must be a segment");
131 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
132 "Dimension of Geometry1 must be higher than that of Geometry2");
133
134 const auto a = geo2.corner(0);
135 const auto b = geo2.corner(1);
136 const auto d = b - a;
137
138 // The initial interval is the whole segment.
139 // Afterwards we start clipping the interval
140 // at the intersections with the facets of geo1
141 tfirst = 0.0;
142 tlast = 1.0;
143
144 const auto facets = getFacetCornerIndices(geo1);
145 for (const auto& f : facets)
146 {
147 const auto n = computeNormal(f);
148
149 const auto c0 = geo1.corner(f[0]);
150 const ctype denom = n*d;
151 const ctype dist = n*(a-c0);
152
153 // use first edge of the facet to scale eps
154 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
155 const ctype eps = baseEps*edge1.two_norm();
156
157 // if denominator is zero the segment in parallel to the edge.
158 // If the distance is positive there is no intersection
159 using std::abs;
160 if (abs(denom) < eps)
161 {
162 if (dist > eps)
163 return false;
164 }
165 else // not parallel: compute line-line intersection
166 {
167 const ctype t = -dist / denom;
168 // if entering half space cut tfirst if t is larger
169 using std::signbit;
170 if (signbit(denom))
171 {
172 if (t > tfirst)
173 tfirst = t;
174 }
175 // if exiting half space cut tlast if t is smaller
176 else
177 {
178 if (t < tlast)
179 tlast = t;
180 }
181 // there is no intersection of the interval is empty
182 // use unscaled epsilon since t is in local coordinates
183 if (tfirst > tlast - baseEps)
184 return false;
185 }
186 }
187
188 // If we made it until here an intersection exists.
189 return true;
190}
191
192} // end namespace Detail
193
200template
201<class Geometry1, class Geometry2,
202 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
203 int dimworld = Geometry1::coorddimension,
204 int dim1 = Geometry1::mydimension,
205 int dim2 = Geometry2::mydimension>
207{
208 static constexpr int dimWorld = Policy::dimWorld;
209
210public:
211 using ctype = typename Policy::ctype;
212 using Point = typename Policy::Point;
213 using Intersection = typename Policy::Intersection;
214
216 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
217 {
218 static_assert(dimworld == Geometry2::coorddimension, "Can only intersect geometries of same coordinate dimension");
219 DUNE_THROW(Dune::NotImplemented, "Geometry intersection detection with intersection computation not implemented for dimworld = "
220 << dimworld << ", dim1 = " << dim1 << ", dim2 = " << dim2);
221 }
222};
223
228template <class Geometry1, class Geometry2, class Policy>
229class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 1>
230{
231 enum { dimworld = 2 };
232 enum { dim1 = 1 };
233 enum { dim2 = 1 };
234
235 // base epsilon for floating point comparisons
236 static constexpr typename Policy::ctype eps_ = 1.5e-7;
237
238public:
239 using ctype = typename Policy::ctype;
240 using Point = typename Policy::Point;
241 using Intersection = typename Policy::Intersection;
242
249 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
250 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
251 {
252 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
253
254 const auto v1 = geo1.corner(1) - geo1.corner(0);
255 const auto v2 = geo2.corner(1) - geo2.corner(0);
256 const auto ac = geo2.corner(0) - geo1.corner(0);
257
258 auto n2 = Point({-1.0*v2[1], v2[0]});
259 n2 /= n2.two_norm();
260
261 // compute distance of first corner on geo2 to line1
262 const auto dist12 = n2*ac;
263
264 // first check parallel segments
265 using std::abs;
266 using std::sqrt;
267
268 const auto v1Norm2 = v1.two_norm2();
269 const auto eps = eps_*sqrt(v1Norm2);
270 const auto eps2 = eps_*v1Norm2;
271
272 const auto sp = n2*v1;
273 if (abs(sp) < eps)
274 {
275 if (abs(dist12) > eps)
276 return false;
277
278 // point intersection can only be one of the corners
279 if ( (v1*v2) < 0.0 )
280 {
281 if ( ac.two_norm2() < eps2 )
282 { intersection = geo2.corner(0); return true; }
283
284 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
285 { intersection = geo2.corner(1); return true; }
286 }
287 else
288 {
289 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
290 { intersection = geo2.corner(1); return true; }
291
292 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
293 { intersection = geo2.corner(0); return true; }
294 }
295
296 // no intersection
297 return false;
298 }
299
300 // intersection point on v1 in local coords
301 const auto t1 = dist12 / sp;
302
303 // check if the local coords are valid
304 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
305 return false;
306
307 // compute global coordinates
308 auto isPoint = geo1.global(t1);
309
310 // check if point is in bounding box of v2
311 const auto c = geo2.corner(0);
312 const auto d = geo2.corner(1);
313
314 using std::min; using std::max;
315 std::array<ctype, 4> bBox({ min(c[0], d[0]), min(c[1], d[1]), max(c[0], d[0]), max(c[1], d[1]) });
316
317 if ( intersectsPointBoundingBox(isPoint, bBox.data()) )
318 {
319 intersection = std::move(isPoint);
320 return true;
321 }
322
323 return false;
324 }
325
332 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
333 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
334 {
335 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
336
337 const auto& a = geo1.corner(0);
338 const auto& b = geo1.corner(1);
339 const auto ab = b - a;
340
341 const auto& p = geo2.corner(0);
342 const auto& q = geo2.corner(1);
343 const auto pq = q - p;
344 Dune::FieldVector<ctype, dimworld> n2{-pq[1], pq[0]};
345
346 using std::max;
347 const auto abNorm2 = ab.two_norm2();
348 const auto pqNorm2 = pq.two_norm2();
349 const auto eps2 = eps_*max(abNorm2, pqNorm2);
350
351 // non-parallel segments do not intersect in a segment.
352 using std::abs;
353 if (abs(n2*ab) > eps2)
354 return false;
355
356 // check if the segments lie on the same line
357 const auto ap = p - a;
358 if (abs(ap*n2) > eps2)
359 return false;
360
361 // compute scaled local coordinates of corner 1/2 of segment2 on segment1
362 auto t1 = ab*ap;
363 auto t2 = ab*(q - a);
364
365 using std::swap;
366 if (t1 > t2)
367 swap(t1, t2);
368
369 using std::clamp;
370 t1 = clamp(t1, 0.0, abNorm2);
371 t2 = clamp(t2, 0.0, abNorm2);
372
373 if (abs(t2-t1) < eps2)
374 return false;
375
376 intersection = Intersection({geo1.global(t1/abNorm2), geo1.global(t2/abNorm2)});
377 return true;
378 }
379};
380
385template <class Geometry1, class Geometry2, class Policy>
386class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 1>
387{
388 enum { dimworld = 2 };
389 enum { dim1 = 2 };
390 enum { dim2 = 1 };
391
392public:
393 using ctype = typename Policy::ctype;
394 using Point = typename Policy::Point;
395 using Intersection = typename Policy::Intersection;
396
397private:
398 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
399
400public:
408 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
409 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
410 {
411 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
412
413 ctype tfirst, tlast;
414 if (intersect_(geo1, geo2, tfirst, tlast))
415 {
416 // the intersection exists. Export the intersections geometry now:
417 // s(t) = a + t(b-a) in [tfirst, tlast]
418 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
419 return true;
420 }
421
422 return false;
423 }
424
432 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
433 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename P::Intersection& intersection)
434 {
435 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
436
437 ctype tfirst, tlast;
438 if (!intersect_(geo1, geo2, tfirst, tlast))
439 {
440 // if start & end point are same, the intersection is a point
441 using std::abs;
442 if ( tfirst > tlast - eps_ )
443 {
444 intersection = Intersection(geo2.global(tfirst));
445 return true;
446 }
447 }
448
449 return false;
450 }
451
452private:
456 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
457 {
458 // lambda to obtain the facet corners on geo1
459 auto getFacetCorners = [] (const Geometry1& geo1)
460 {
461 std::vector< std::array<int, 2> > facetCorners;
462 switch (geo1.corners())
463 {
464 case 4: // quadrilateral
465 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
466 break;
467 case 3: // triangle
468 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
469 break;
470 default:
471 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
472 << geo1.type() << " with "<< geo1.corners() << " corners.");
473 }
474
475 return facetCorners;
476 };
477
478 // lambda to obtain the normal on a facet
479 const auto center1 = geo1.center();
480 auto computeNormal = [&geo1, &center1] (const std::array<int, 2>& facetCorners)
481 {
482 const auto c0 = geo1.corner(facetCorners[0]);
483 const auto c1 = geo1.corner(facetCorners[1]);
484 const auto edge = c1 - c0;
485
486 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
487 n /= n.two_norm();
488
489 // orientation of the element is not known
490 // make sure normal points outwards of element
491 if ( n*(center1-c0) > 0.0 )
492 n *= -1.0;
493
494 return n;
495 };
496
497 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
498 }
499};
500
505template <class Geometry1, class Geometry2, class Policy>
506class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 2>
507: public GeometryIntersection<Geometry2, Geometry1, Policy, 2, 2, 1>
508{
510
511public:
519 template<class P = Policy>
520 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
521 {
522 return Base::intersection(geo2, geo1, intersection);
523 }
524};
525
530template <class Geometry1, class Geometry2, class Policy>
531class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 2>
532{
533 enum { dimworld = 2 };
534 enum { dim1 = 2 };
535 enum { dim2 = 2 };
536
537public:
538 using ctype = typename Policy::ctype;
539 using Point = typename Policy::Point;
540 using Intersection = typename Policy::Intersection;
541
542private:
543 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
544
545public:
558 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
559 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
560 {
561 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
562
563 // the candidate intersection points
564 std::vector<Point> points; points.reserve(6);
565
566 // add polygon1 corners that are inside polygon2
567 for (int i = 0; i < geo1.corners(); ++i)
568 if (intersectsPointGeometry(geo1.corner(i), geo2))
569 points.emplace_back(geo1.corner(i));
570
571 const auto numPoints1 = points.size();
572 if (numPoints1 != geo1.corners())
573 {
574 // add polygon2 corners that are inside polygon1
575 for (int i = 0; i < geo2.corners(); ++i)
576 if (intersectsPointGeometry(geo2.corner(i), geo1))
577 points.emplace_back(geo2.corner(i));
578
579 if (points.empty())
580 return false;
581
582 if (points.size() - numPoints1 != geo2.corners())
583 {
584 const auto refElement1 = referenceElement(geo1);
585 const auto refElement2 = referenceElement(geo2);
586
587 // add intersections of edges
588 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
590 for (int i = 0; i < refElement1.size(dim1-1); ++i)
591 {
592 const auto localEdgeGeom1 = refElement1.template geometry<dim1-1>(i);
593 const auto edge1 = SegGeometry( Dune::GeometryTypes::line,
594 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
595 geo1.global(localEdgeGeom1.corner(1))} ));
596
597 for (int j = 0; j < refElement2.size(dim2-1); ++j)
598 {
599 const auto localEdgeGeom2 = refElement2.template geometry<dim2-1>(j);
600 const auto edge2 = SegGeometry( Dune::GeometryTypes::line,
601 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
602 geo2.global(localEdgeGeom2.corner(1))} ));
603
605 typename EdgeTest::Intersection edgeIntersection;
606 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
607 points.emplace_back(edgeIntersection);
608 }
609 }
610 }
611 }
612
613 if (points.empty())
614 return false;
615
616 // remove duplicates
617 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
618 std::sort(points.begin(), points.end(), [&eps](const auto& a, const auto& b) -> bool
619 {
620 using std::abs;
621 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
622 });
623
624 auto removeIt = std::unique(points.begin(), points.end(), [&eps](const auto& a, const auto&b)
625 {
626 return (b-a).two_norm() < eps;
627 });
628
629 points.erase(removeIt, points.end());
630
631 // return false if we don't have at least three unique points
632 if (points.size() < 3)
633 return false;
634
635 // intersection polygon is convex hull of above points
636 intersection = grahamConvexHull<2>(points);
637 assert(!intersection.empty());
638 return true;
639 }
640
648 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
649 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
650 {
651 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
652 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for segment-like intersections");
653 }
654
662 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
663 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
664 {
665 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
666 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for touching points");
667 }
668};
669
674template <class Geometry1, class Geometry2, class Policy>
675class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 1>
676{
677 enum { dimworld = 3 };
678 enum { dim1 = 3 };
679 enum { dim2 = 1 };
680
681public:
682 using ctype = typename Policy::ctype;
683 using Point = typename Policy::Point;
684 using Intersection = typename Policy::Intersection;
685
686private:
687 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
688
689public:
701 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
702 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
703 {
704 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
705
706 ctype tfirst, tlast;
707 if (intersect_(geo1, geo2, tfirst, tlast))
708 {
709 // Intersection exists. Export the intersections geometry now:
710 // s(t) = a + t(b-a) in [tfirst, tlast]
711 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
712 return true;
713 }
714
715 return false;
716 }
717
729 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
730 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
731 {
732 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
733
734 ctype tfirst, tlast;
735 if (intersect_(geo1, geo2, tfirst, tlast))
736 {
737 // if start & end point are same, the intersection is a point
738 using std::abs;
739 if ( abs(tfirst - tlast) < eps_ )
740 {
741 intersection = Intersection(geo2.global(tfirst));
742 return true;
743 }
744 }
745
746 return false;
747 }
748
749private:
753 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
754 {
755 // lambda to obtain facet corners on geo1
756 auto getFacetCorners = [] (const Geometry1& geo1)
757 {
758 std::vector< std::vector<int> > facetCorners;
759 // sort facet corner so that normal n = (p1-p0)x(p2-p0) always points outwards
760 switch (geo1.corners())
761 {
762 // todo compute intersection geometries
763 case 8: // hexahedron
764 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
765 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
766 break;
767 case 6: // prism
768 facetCorners = {{0, 2, 1}, {3, 4, 5}, {0, 1, 3, 4}, {0, 3, 2, 5}, {1, 2, 4, 5}};
769 break;
770 case 5: // pyramid
771 facetCorners = {{0, 2, 1, 3}, {0, 1, 4}, {1, 3, 4}, {2, 4, 3}, {0, 4, 2}};
772 break;
773 case 4: // tetrahedron
774 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
775 break;
776 default:
777 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
778 << geo1.type() << ", "<< geo1.corners() << " corners.");
779 }
780
781 return facetCorners;
782 };
783
784 // lambda to obtain the normal on a facet
785 auto computeNormal = [&geo1] (const std::vector<int>& facetCorners)
786 {
787 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
788 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
789
790 auto n = crossProduct(v0, v1);
791 n /= n.two_norm();
792
793 return n;
794 };
795
796 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
797 }
798};
799
804template <class Geometry1, class Geometry2, class Policy>
805class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 3>
806: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 1>
807{
809public:
817 template<class P = Policy>
818 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
819 {
820 return Base::intersection(geo2, geo1, intersection);
821 }
822};
823
828template <class Geometry1, class Geometry2, class Policy>
829class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 2>
830{
831 enum { dimworld = 3 };
832 enum { dim1 = 2 };
833 enum { dim2 = 2 };
834
835public:
836 using ctype = typename Policy::ctype;
837 using Point = typename Policy::Point;
838 using Intersection = typename Policy::Intersection;
839
840private:
841 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
842
843public:
855 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
856 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
857 {
858 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
859
860 // if the geometries are not parallel, there cannot be a polygon intersection
861 const auto a1 = geo1.corner(1) - geo1.corner(0);
862 const auto b1 = geo1.corner(2) - geo1.corner(0);
863 const auto n1 = crossProduct(a1, b1);
864
865 const auto a2 = geo2.corner(1) - geo2.corner(0);
866 const auto b2 = geo2.corner(2) - geo2.corner(0);
867 const auto n2 = crossProduct(a2, b2);
868
869 // compute an epsilon scaled with larger edge length
870 const auto a1Norm2 = a1.two_norm2();
871 const auto b1Norm2 = b1.two_norm2();
872
873 using std::max;
874 const auto maxNorm2 = max(a1Norm2, b1Norm2);
875 const auto eps2 = maxNorm2*eps_;
876
877 // early exit if polygons are not parallel
878 using std::abs;
879 if (crossProduct(n1, n2).two_norm2() > eps2*maxNorm2)
880 return false;
881
882 // they are parallel, verify that they are on the same plane
883 auto d1 = geo2.corner(0) - geo1.corner(0);
884 if (d1.two_norm2() < eps2)
885 d1 = geo2.corner(1) - geo1.corner(0);
886
887 using std::sqrt;
888 const auto maxNorm = sqrt(maxNorm2);
889 const auto eps3 = maxNorm*eps2;
890
891 if (abs(d1*n2) > eps3)
892 return false;
893
894 // the candidate intersection points
895 std::vector<Point> points; points.reserve(8);
896
897 // add polygon1 corners that are inside polygon2
898 for (int i = 0; i < geo1.corners(); ++i)
899 if (intersectsPointGeometry(geo1.corner(i), geo2))
900 points.emplace_back(geo1.corner(i));
901
902 const auto numPoints1 = points.size();
903 const bool resultIsGeo1 = numPoints1 == geo1.corners();
904 if (!resultIsGeo1)
905 {
906 // add polygon2 corners that are inside polygon1
907 for (int i = 0; i < geo2.corners(); ++i)
908 if (intersectsPointGeometry(geo2.corner(i), geo1))
909 points.emplace_back(geo2.corner(i));
910
911 const bool resultIsGeo2 = (points.size() - numPoints1) == geo2.corners();
912 if (!resultIsGeo2)
913 {
914 const auto referenceElement1 = referenceElement(geo1);
915 const auto referenceElement2 = referenceElement(geo2);
916
917 // add intersections of edges
918 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
920 for (int i = 0; i < referenceElement1.size(dim1-1); ++i)
921 {
922 const auto localEdgeGeom1 = referenceElement1.template geometry<dim1-1>(i);
923 const auto edge1 = SegGeometry(
925 std::vector<Point>({
926 geo1.global(localEdgeGeom1.corner(0)),
927 geo1.global(localEdgeGeom1.corner(1))
928 })
929 );
930
931 for (int j = 0; j < referenceElement2.size(dim2-1); ++j)
932 {
933 const auto localEdgeGeom2 = referenceElement2.template geometry<dim2-1>(j);
934 const auto edge2 = SegGeometry(
936 std::vector<Point>({
937 geo2.global(localEdgeGeom2.corner(0)),
938 geo2.global(localEdgeGeom2.corner(1))
939 })
940 );
941
943 typename EdgeTest::Intersection edgeIntersection;
944 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
945 points.emplace_back(edgeIntersection);
946 }
947 }
948 }
949 }
950
951 if (points.empty())
952 return false;
953
954 // remove duplicates
955 const auto eps = maxNorm*eps_;
956 std::sort(points.begin(), points.end(), [eps] (const auto& a, const auto& b) -> bool
957 {
958 using std::abs;
959 return (abs(a[0]-b[0]) > eps ? a[0] < b[0]
960 : (abs(a[1]-b[1]) > eps ? a[1] < b[1]
961 : (a[2] < b[2])));
962 });
963
964 const auto squaredEps = eps*eps;
965 points.erase(std::unique(
966 points.begin(), points.end(),
967 [squaredEps] (const auto& a, const auto&b) { return (b-a).two_norm2() < squaredEps; }),
968 points.end()
969 );
970
971 // return false if we don't have at least three unique points
972 if (points.size() < 3)
973 return false;
974
975 // intersection polygon is convex hull of above points
976 intersection = grahamConvexHull<2>(points);
977 assert(!intersection.empty());
978 return true;
979 }
980
988 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
989 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
990 {
991 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
992 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for segment-like intersections");
993 }
994
1002 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1003 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1004 {
1005 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1006 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for touching points");
1007 }
1008};
1009
1014template <class Geometry1, class Geometry2, class Policy>
1015class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 2>
1016{
1017 enum { dimworld = 3 };
1018 enum { dim1 = 3 };
1019 enum { dim2 = 2 };
1020
1021public:
1022 using ctype = typename Policy::ctype;
1023 using Point = typename Policy::Point;
1024 using Intersection = typename Policy::Intersection;
1025
1026private:
1027 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1028
1029public:
1044 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
1045 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1046 {
1047 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1048
1049 // the candidate intersection points
1050 std::vector<Point> points; points.reserve(10);
1051
1052 // add 3d geometry corners that are inside the 2d geometry
1053 for (int i = 0; i < geo1.corners(); ++i)
1054 if (intersectsPointGeometry(geo1.corner(i), geo2))
1055 points.emplace_back(geo1.corner(i));
1056
1057 // add 2d geometry corners that are inside the 3d geometry
1058 for (int i = 0; i < geo2.corners(); ++i)
1059 if (intersectsPointGeometry(geo2.corner(i), geo1))
1060 points.emplace_back(geo2.corner(i));
1061
1062 // get some geometry types
1063 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
1064 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
1065
1066 const auto refElement1 = referenceElement(geo1);
1067 const auto refElement2 = referenceElement(geo2);
1068
1069 // intersection policy for point-like intersections (used later)
1071
1072 // add intersection points of all polyhedron edges (codim dim-1) with the polygon
1073 for (int i = 0; i < refElement1.size(dim1-1); ++i)
1074 {
1075 const auto localEdgeGeom = refElement1.template geometry<dim1-1>(i);
1076 const auto p = geo1.global(localEdgeGeom.corner(0));
1077 const auto q = geo1.global(localEdgeGeom.corner(1));
1078 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
1079
1081 typename PolySegTest::Intersection polySegIntersection;
1082 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
1083 points.emplace_back(std::move(polySegIntersection));
1084 }
1085
1086 // add intersection points of all polygon faces (codim 1) with the polyhedron faces
1087 for (int i = 0; i < refElement1.size(1); ++i)
1088 {
1089 const auto faceGeo = [&]()
1090 {
1091 const auto localFaceGeo = refElement1.template geometry<1>(i);
1092 if (localFaceGeo.corners() == 4)
1093 {
1094 const auto a = geo1.global(localFaceGeo.corner(0));
1095 const auto b = geo1.global(localFaceGeo.corner(1));
1096 const auto c = geo1.global(localFaceGeo.corner(2));
1097 const auto d = geo1.global(localFaceGeo.corner(3));
1098
1099 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
1100 }
1101 else
1102 {
1103 const auto a = geo1.global(localFaceGeo.corner(0));
1104 const auto b = geo1.global(localFaceGeo.corner(1));
1105 const auto c = geo1.global(localFaceGeo.corner(2));
1106
1107 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
1108 }
1109 }();
1110
1111 for (int j = 0; j < refElement2.size(1); ++j)
1112 {
1113 const auto localEdgeGeom = refElement2.template geometry<1>(j);
1114 const auto p = geo2.global(localEdgeGeom.corner(0));
1115 const auto q = geo2.global(localEdgeGeom.corner(1));
1116
1117 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
1118
1120 typename PolySegTest::Intersection polySegIntersection;
1121 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
1122 points.emplace_back(std::move(polySegIntersection));
1123 }
1124 }
1125
1126 // return if no intersection points were found
1127 if (points.empty()) return false;
1128
1129 // remove duplicates
1130 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
1131 const auto notEqual = [eps] (auto a, auto b) { using std::abs; return abs(b-a) > eps; };
1132 std::sort(points.begin(), points.end(), [notEqual](const auto& a, const auto& b) -> bool
1133 {
1134 return (notEqual(a[0], b[0]) ? a[0] < b[0]
1135 : (notEqual(a[1], b[1]) ? a[1] < b[1]
1136 : (a[2] < b[2])));
1137 });
1138
1139 const auto squaredEps = eps*eps;
1140 points.erase(std::unique(
1141 points.begin(), points.end(),
1142 [squaredEps] (const auto& a, const auto&b) { return (b-a).two_norm2() < squaredEps; }),
1143 points.end()
1144 );
1145
1146 // return false if we don't have more than three unique points
1147 if (points.size() < 3) return false;
1148
1149 // intersection polygon is convex hull of above points
1150 intersection = grahamConvexHull<2>(points);
1151 assert(!intersection.empty());
1152
1153 return true;
1154 }
1155
1170 template<class P = Policy, std::enable_if_t<P::dimIntersection != 2, int> = 0>
1171 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1172 {
1173 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1174 DUNE_THROW(Dune::NotImplemented, "Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
1175 }
1176};
1177
1182template <class Geometry1, class Geometry2, class Policy>
1183class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 3>
1184: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 2>
1185{
1187public:
1195 template<class P = Policy>
1196 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
1197 {
1198 return Base::intersection(geo2, geo1, intersection);
1199 }
1200};
1201
1206template <class Geometry1, class Geometry2, class Policy>
1207class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 3>
1208{
1209 enum { dimworld = 3 };
1210 enum { dim1 = 3 };
1211 enum { dim2 = 3 };
1212
1213public:
1214 using ctype = typename Policy::ctype;
1215 using Point = typename Policy::Point;
1216 using Intersection = typename Policy::Intersection;
1217
1218private:
1219 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1220
1221public:
1233 template<class P = Policy, std::enable_if_t<P::dimIntersection == 3, int> = 0>
1234 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1235 {
1236 static_assert(
1237 int(dimworld) == int(Geometry2::coorddimension),
1238 "Can only collide geometries of same coordinate dimension"
1239 );
1240
1241 const auto refElement1 = referenceElement(geo1);
1242 const auto refElement2 = referenceElement(geo2);
1243
1244 // the candidate intersection points
1245 std::vector<Point> points;
1246 points.reserve(refElement1.size(2) + refElement2.size(2));
1247
1248 // add corners inside the other geometry
1249 const auto addPointIntersections = [&](const auto& g1, const auto& g2)
1250 {
1251 for (int i = 0; i < g1.corners(); ++i)
1252 if (const auto& c = g1.corner(i); intersectsPointGeometry(c, g2))
1253 points.emplace_back(c);
1254 };
1255
1256 addPointIntersections(geo1, geo2);
1257 addPointIntersections(geo2, geo1);
1258
1259 // get geometry types for the facets
1260 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
1261 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
1262
1263 // intersection policy for point-like intersections (used later)
1265
1266 // add intersection points of all edges with the faces of the other polyhedron
1267 const auto addEdgeIntersections = [&](const auto& g1, const auto& g2, const auto& ref1, const auto& ref2)
1268 {
1269 for (int i = 0; i < ref1.size(1); ++i)
1270 {
1271 const auto faceGeo = [&]()
1272 {
1273 const auto localFaceGeo = ref1.template geometry<1>(i);
1274 if (localFaceGeo.corners() == 4)
1275 {
1276 const auto a = g1.global(localFaceGeo.corner(0));
1277 const auto b = g1.global(localFaceGeo.corner(1));
1278 const auto c = g1.global(localFaceGeo.corner(2));
1279 const auto d = g1.global(localFaceGeo.corner(3));
1280
1281 return PolyhedronFaceGeometry(
1282 Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d}
1283 );
1284 }
1285 else
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
1291 return PolyhedronFaceGeometry(
1292 Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c}
1293 );
1294 }
1295 }();
1296
1297 for (int j = 0; j < ref2.size(2); ++j)
1298 {
1299 const auto localEdgeGeom = ref2.template geometry<2>(j);
1300 const auto p = g2.global(localEdgeGeom.corner(0));
1301 const auto q = g2.global(localEdgeGeom.corner(1));
1302
1303 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
1304
1306 typename PolySegTest::Intersection polySegIntersection;
1307 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
1308 points.emplace_back(std::move(polySegIntersection));
1309 }
1310 }
1311 };
1312
1313 addEdgeIntersections(geo1, geo2, refElement1, refElement2);
1314 addEdgeIntersections(geo2, geo1, refElement2, refElement1);
1315
1316 // return if no intersection points were found
1317 if (points.empty())
1318 return false;
1319
1320 // remove duplicates
1321 const auto norm = (geo1.corner(0) - geo1.corner(1)).two_norm();
1322 const auto eps = norm*eps_;
1323 const auto notEqual = [eps] (auto a, auto b) { using std::abs; return abs(b-a) > eps; };
1324 std::sort(points.begin(), points.end(), [notEqual](const auto& a, const auto& b) -> bool
1325 {
1326 return (notEqual(a[0], b[0]) ? a[0] < b[0]
1327 : (notEqual(a[1], b[1]) ? a[1] < b[1]
1328 : (a[2] < b[2])));
1329 });
1330
1331 const auto squaredEps = eps*eps;
1332 points.erase(std::unique(
1333 points.begin(), points.end(),
1334 [squaredEps] (const auto& a, const auto&b) { return (b-a).two_norm2() < squaredEps; }),
1335 points.end()
1336 );
1337
1338 // return false if we don't have more than four unique points (dim+1)
1339 if (points.size() < 4)
1340 return false;
1341
1342 // check if the points are coplanar (then we don't have a 3-dimensional intersection)
1343 const bool coplanar = [&]
1344 {
1345 const auto p0 = points[0];
1346 const auto normal = crossProduct(points[1] - p0, points[2] - p0);
1347 // include the normal in eps (instead of norm*norm) since the normal can become very small
1348 // if the first three points are very close together
1349 const auto epsCoplanar = normal.two_norm()*norm*eps_;
1350 for (int i = 3; i < points.size(); ++i)
1351 {
1352 const auto ad = points[i] - p0;
1353 using std::abs;
1354 if (abs(normal*ad) > epsCoplanar)
1355 return false;
1356 }
1357
1358 return true;
1359 }();
1360
1361 if (coplanar)
1362 return false;
1363
1364 // we return the intersection as point cloud (that can be triangulated later)
1365 intersection = points;
1366
1367 return true;
1368 }
1369
1376 template<class P = Policy, std::enable_if_t<P::dimIntersection != 3, int> = 0>
1377 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1378 {
1379 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1380 DUNE_THROW(Dune::NotImplemented, "Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
1381 }
1382};
1383
1388template <class Geometry1, class Geometry2, class Policy>
1389class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 1>
1390{
1391 enum { dimworld = 3 };
1392 enum { dim1 = 2 };
1393 enum { dim2 = 1 };
1394
1395public:
1396 using ctype = typename Policy::ctype;
1397 using Point = typename Policy::Point;
1398 using Intersection = typename Policy::Intersection;
1399
1400private:
1401 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1402
1403public:
1413 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1414 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1415 {
1416 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1417
1418 ctype tfirst, tlast;
1419 if (intersect_(geo1, geo2, tfirst, tlast))
1420 {
1421 // the intersection exists. Export the intersections geometry now:
1422 // s(t) = a + t(b-a) in [tfirst, tlast]
1423 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
1424 return true;
1425 }
1426
1427 return false;
1428 }
1429
1439 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1440 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& is)
1441 {
1442 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1443
1444 const auto p = geo2.corner(0);
1445 const auto q = geo2.corner(1);
1446
1447 const auto a = geo1.corner(0);
1448 const auto b = geo1.corner(1);
1449 const auto c = geo1.corner(2);
1450
1451 if (geo1.corners() == 3)
1452 return intersect_<Policy>(a, b, c, p, q, is);
1453
1454 else if (geo1.corners() == 4)
1455 {
1456 Intersection is1, is2;
1457 bool hasSegment1, hasSegment2;
1458
1459 const auto d = geo1.corner(3);
1460 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1461 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1462
1463 if (hasSegment1 || hasSegment2)
1464 return false;
1465
1466 if (intersects1) { is = is1; return true; }
1467 if (intersects2) { is = is2; return true; }
1468
1469 return false;
1470 }
1471
1472 else
1473 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1474 << geo1.type() << ", "<< geo1.corners() << " corners.");
1475 }
1476
1477private:
1481 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1482 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
1483 {
1484 // lambda to obtain the facet corners on geo1
1485 auto getFacetCorners = [] (const Geometry1& geo1)
1486 {
1487 std::vector< std::array<int, 2> > facetCorners;
1488 switch (geo1.corners())
1489 {
1490 case 4: // quadrilateral
1491 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1492 break;
1493 case 3: // triangle
1494 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1495 break;
1496 default:
1497 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1498 << geo1.type() << " with "<< geo1.corners() << " corners.");
1499 }
1500
1501 return facetCorners;
1502 };
1503
1504 const auto center1 = geo1.center();
1505 const auto normal1 = crossProduct(geo1.corner(1) - geo1.corner(0),
1506 geo1.corner(2) - geo1.corner(0));
1507
1508 // lambda to obtain the normal on a facet
1509 auto computeNormal = [&center1, &normal1, &geo1] (const std::array<int, 2>& facetCorners)
1510 {
1511 const auto c0 = geo1.corner(facetCorners[0]);
1512 const auto c1 = geo1.corner(facetCorners[1]);
1513 const auto edge = c1 - c0;
1514
1515 Dune::FieldVector<ctype, dimworld> n = crossProduct(edge, normal1);
1516 n /= n.two_norm();
1517
1518 // orientation of the element is not known
1519 // make sure normal points outwards of element
1520 if ( n*(center1-c0) > 0.0 )
1521 n *= -1.0;
1522
1523 return n;
1524 };
1525
1526 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
1527 }
1528
1532 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1533 static bool intersect_(const Point& a, const Point& b, const Point& c,
1534 const Point& p, const Point& q,
1535 Intersection& is)
1536 {
1537 bool hasSegment;
1538 return intersect_(a, b, c, p, q, is, hasSegment);
1539 }
1540
1545 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1546 static bool intersect_(const Point& a, const Point& b, const Point& c,
1547 const Point& p, const Point& q,
1548 Intersection& is, bool& hasSegmentIs)
1549 {
1550 hasSegmentIs = false;
1551
1552 const auto ab = b - a;
1553 const auto ac = c - a;
1554 const auto qp = p - q;
1555
1556 // compute the triangle normal that defines the triangle plane
1557 const auto normal = crossProduct(ab, ac);
1558
1559 // compute the denominator
1560 // if denom is 0 the segment is parallel and we can return
1561 const auto denom = normal*qp;
1562 const auto abnorm2 = ab.two_norm2();
1563 const auto eps = eps_*abnorm2*qp.two_norm();
1564 using std::abs;
1565 if (abs(denom) < eps)
1566 {
1567 const auto pa = a - p;
1568 const auto denom2 = normal*pa;
1569 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1570 return false;
1571
1572 // if we get here, we are in-plane. Check if a
1573 // segment-like intersection with geometry 1 exists.
1574 using SegmentPolicy = typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1575 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1576 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1577 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1578 using SegmentIntersectionType = typename SegmentIntersectionAlgorithm::Intersection;
1579 SegmentIntersectionType segmentIs;
1580
1581 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1582 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1583 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1584 {
1585 hasSegmentIs = true;
1586 return false;
1587 }
1588
1589 // We are in-plane. A point-like
1590 // intersection can only be on the edges
1591 for (const auto& ip : {p, q})
1592 {
1593 if ( intersectsPointSimplex(ip, a, b)
1594 || intersectsPointSimplex(ip, b, c)
1595 || intersectsPointSimplex(ip, c, a) )
1596 {
1597 is = ip;
1598 return true;
1599 }
1600 }
1601
1602 return false;
1603 }
1604
1605 // compute intersection t value of pq with plane of triangle.
1606 // a segment intersects if and only if 0 <= t <= 1.
1607 const auto ap = p - a;
1608 const auto t = (ap*normal)/denom;
1609 if (t < 0.0 - eps_) return false;
1610 if (t > 1.0 + eps_) return false;
1611
1612 // compute the barycentric coordinates and check if the intersection point
1613 // is inside the bounds of the triangle
1614 const auto e = crossProduct(qp, ap);
1615 const auto v = (ac*e)/denom;
1616 if (v < -eps_ || v > 1.0 + eps_) return false;
1617 const auto w = -(ab*e)/denom;
1618 if (w < -eps_ || v + w > 1.0 + eps_) return false;
1619
1620 // Now we are sure there is an intersection points
1621 // Perform delayed division compute the last barycentric coordinate component
1622 const auto u = 1.0 - v - w;
1623
1624 Point ip(0.0);
1625 ip.axpy(u, a);
1626 ip.axpy(v, b);
1627 ip.axpy(w, c);
1628 is = ip;
1629 return true;
1630 }
1631};
1632
1637template <class Geometry1, class Geometry2, class Policy>
1638class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 2>
1639: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 2, 1>
1640{
1642public:
1650 template<class P = Policy>
1651 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
1652 {
1653 return Base::intersection(geo2, geo1, intersection);
1654 }
1655};
1656
1661template <class Geometry1, class Geometry2, class Policy>
1662class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 1>
1663{
1664 enum { dimworld = 3 };
1665 enum { dim1 = 1 };
1666 enum { dim2 = 1 };
1667
1668public:
1669 using ctype = typename Policy::ctype;
1670 using Point = typename Policy::Point;
1671 using Intersection = typename Policy::Intersection;
1672
1673private:
1674 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1675
1676public:
1683 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1684 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1685 {
1686 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1687
1688 const auto v1 = geo1.corner(1) - geo1.corner(0);
1689 const auto v2 = geo2.corner(1) - geo2.corner(0);
1690 const auto ac = geo2.corner(0) - geo1.corner(0);
1691
1692 const auto v1Norm2 = v1.two_norm2();
1693 const auto eps2 = eps_*v1Norm2;
1694
1695 const auto n = crossProduct(v1, v2);
1696
1697 // first check if segments are parallel
1698 using std::abs;
1699 if ( n.two_norm2() < eps2*v1Norm2 )
1700 {
1701 // check if they lie on the same line
1702 if (crossProduct(v1, ac).two_norm2() > eps2)
1703 return false;
1704
1705 // they lie on the same line,
1706 // if so, point intersection can only be one of the corners
1707 const auto sp = v1*v2;
1708 if ( sp < 0.0 )
1709 {
1710 if ( ac.two_norm2() < eps2 )
1711 { intersection = geo2.corner(0); return true; }
1712
1713 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
1714 { intersection = geo2.corner(1); return true; }
1715 }
1716 else
1717 {
1718 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
1719 { intersection = geo2.corner(1); return true; }
1720
1721 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
1722 { intersection = geo2.corner(0); return true; }
1723 }
1724
1725 // no intersection
1726 return false;
1727 }
1728
1729 // in-plane normal vector on v1
1730 const auto v1Normal = crossProduct(v1, n);
1731
1732 // intersection point on v2 in local coords
1733 const auto t2 = - 1.0*(ac*v1Normal) / (v2*v1Normal);
1734
1735 // check if the local coords are valid
1736 if (t2 < -1.0*eps_ || t2 > 1.0 + eps_)
1737 return false;
1738
1739 if (auto ip = geo2.global(t2); intersectsPointGeometry(ip, geo1))
1740 { intersection = std::move(ip); return true; }
1741
1742 return false;
1743 }
1744
1752 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1753 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1754 {
1755 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1756
1757 const auto& a = geo1.corner(0);
1758 const auto& b = geo1.corner(1);
1759 const auto ab = b-a;
1760
1761 const auto& p = geo2.corner(0);
1762 const auto& q = geo2.corner(1);
1763 const auto pq = q-p;
1764
1765 const auto abNorm2 = ab.two_norm2();
1766 const auto pqNorm2 = pq.two_norm2();
1767
1768 using std::max;
1769 const auto eps2 = eps_*max(abNorm2, pqNorm2);
1770
1771 // if the segment are not parallel there is no segment intersection
1772 using std::abs;
1773 if (crossProduct(ab, pq).two_norm2() > eps2*eps2)
1774 return false;
1775
1776 const auto ap = (p-a);
1777 const auto aq = (q-a);
1778
1779 // points have to be colinear
1780 if (crossProduct(ap, aq).two_norm2() > eps2*eps2)
1781 return false;
1782
1783 // scaled local coordinates
1784 // we save the division for the normalization for now
1785 // and do it in the very end if we find an intersection
1786 auto tp = ap*ab;
1787 auto tq = aq*ab;
1788
1789 // make sure they are sorted
1790 using std::swap;
1791 if (tp > tq)
1792 swap(tp, tq);
1793
1794 using std::clamp;
1795 tp = clamp(tp, 0.0, abNorm2);
1796 tq = clamp(tq, 0.0, abNorm2);
1797
1798 if (abs(tp-tq) < eps2)
1799 return false;
1800
1801 intersection = Intersection({geo1.global(tp/abNorm2), geo1.global(tq/abNorm2)});
1802 return true;
1803 }
1804};
1805
1806} // end namespace Dumux
1807
1808#endif
An axis-aligned bounding box volume hierarchy for dune grids.
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:241
typename Policy::ctype ctype
Definition: geometryintersection.hh:239
typename Policy::Point Point
Definition: geometryintersection.hh:240
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:250
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:520
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:409
typename Policy::ctype ctype
Definition: geometryintersection.hh:393
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:395
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:433
typename Policy::Point Point
Definition: geometryintersection.hh:394
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:540
typename Policy::ctype ctype
Definition: geometryintersection.hh:538
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:559
typename Policy::Point Point
Definition: geometryintersection.hh:539
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1671
typename Policy::Point Point
Definition: geometryintersection.hh:1670
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:1684
typename Policy::ctype ctype
Definition: geometryintersection.hh:1669
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:1651
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:818
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1398
typename Policy::Point Point
Definition: geometryintersection.hh:1397
typename Policy::ctype ctype
Definition: geometryintersection.hh:1396
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1414
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1440
typename Policy::Point Point
Definition: geometryintersection.hh:837
typename Policy::ctype ctype
Definition: geometryintersection.hh:836
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:838
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:856
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:1196
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:684
typename Policy::Point Point
Definition: geometryintersection.hh:683
typename Policy::ctype ctype
Definition: geometryintersection.hh:682
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:702
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1024
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:1045
typename Policy::Point Point
Definition: geometryintersection.hh:1023
typename Policy::ctype ctype
Definition: geometryintersection.hh:1022
typename Policy::Point Point
Definition: geometryintersection.hh:1215
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1216
typename Policy::ctype ctype
Definition: geometryintersection.hh:1214
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two convex polyhedra.
Definition: geometryintersection.hh:1234
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:207
typename Policy::ctype ctype
Definition: geometryintersection.hh:211
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:213
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Determine if the two geometries intersect and compute the intersection geometry.
Definition: geometryintersection.hh:216
typename Policy::Point Point
Definition: geometryintersection.hh:212
default policy chooser class
Definition: geometryintersection.hh:85
std::tuple_element_t< isDim, DefaultPolicies > type
Definition: geometryintersection.hh:99
A function to compute the convex hull of a point cloud.
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:671
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:28
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: normal.hh:26
bool intersectsPointSimplex(const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &p0, const Dune::FieldVector< ctype, dimworld > &p1, const Dune::FieldVector< ctype, dimworld > &p2, const Dune::FieldVector< ctype, dimworld > &p3)
Find out whether a point is inside the tetrahedron (p0, p1, p2, p3) (dimworld is 3)
Definition: intersectspointsimplex.hh:26
bool computeSegmentIntersection(const Geo1 &geo1, const Geo2 &geo2, ctype baseEps, ctype &tfirst, ctype &tlast, const GetFacetCornerIndices &getFacetCornerIndices, const ComputeNormalFunction &computeNormal)
Algorithm to find segment-like intersections of a polygon/polyhedron with a segment....
Definition: geometryintersection.hh:125
Detect if a point intersects a geometry.
Define some often used mathematical functions.
constexpr Line line
Definition: couplingmanager1d3d_line.hh:31
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:104
Definition: adapt.hh:17
bool intersectsPointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Check whether a point is intersectin a bounding box (dimworld == 3)
Definition: boundingboxtree.hh:294
Policy structure for point-like intersections.
Definition: geometryintersection.hh:34
static constexpr int dimIntersection
Definition: geometryintersection.hh:39
Point Intersection
Definition: geometryintersection.hh:40
ct ctype
Definition: geometryintersection.hh:35
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:36
static constexpr int dimWorld
Definition: geometryintersection.hh:38
Policy structure for polygon-like intersections.
Definition: geometryintersection.hh:59
static constexpr int dimWorld
Definition: geometryintersection.hh:64
static constexpr int dimIntersection
Definition: geometryintersection.hh:65
ct ctype
Definition: geometryintersection.hh:60
Polygon Intersection
Definition: geometryintersection.hh:66
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:61
std::vector< Point > Polygon
Definition: geometryintersection.hh:62
Policy structure for polyhedron-like intersections.
Definition: geometryintersection.hh:72
std::vector< Point > Polyhedron
Definition: geometryintersection.hh:75
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:74
Polyhedron Intersection
Definition: geometryintersection.hh:79
static constexpr int dimIntersection
Definition: geometryintersection.hh:78
static constexpr int dimWorld
Definition: geometryintersection.hh:77
ct ctype
Definition: geometryintersection.hh:73
Policy structure for segment-like intersections.
Definition: geometryintersection.hh:46
Segment Intersection
Definition: geometryintersection.hh:53
static constexpr int dimIntersection
Definition: geometryintersection.hh:52
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:48
static constexpr int dimWorld
Definition: geometryintersection.hh:51
std::array< Point, 2 > Segment
Definition: geometryintersection.hh:49
ct ctype
Definition: geometryintersection.hh:47