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