3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
geometry/geometryintersection.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
23#ifndef DUMUX_GEOMETRY_INTERSECTION_HH
24#define DUMUX_GEOMETRY_INTERSECTION_HH
25
26#include <tuple>
27
28#include <dune/common/exceptions.hh>
29#include <dune/common/promotiontraits.hh>
30#include <dune/geometry/multilineargeometry.hh>
31
32#include <dumux/common/math.hh>
36
37namespace Dumux {
38namespace IntersectionPolicy {
39
41template<class ct, int dw>
43{
44 using ctype = ct;
45 using Point = Dune::FieldVector<ctype, dw>;
46
47 static constexpr int dimWorld = dw;
48 static constexpr int dimIntersection = 0;
50};
51
53template<class ct, int dw>
55{
56 using ctype = ct;
57 using Point = Dune::FieldVector<ctype, dw>;
58 using Segment = std::array<Point, 2>;
59
60 static constexpr int dimWorld = dw;
61 static constexpr int dimIntersection = 1;
63};
64
66template<class ct, int dw>
68{
69 using ctype = ct;
70 using Point = Dune::FieldVector<ctype, dw>;
71 using Polygon = std::vector<Point>;
72
73 static constexpr int dimWorld = dw;
74 static constexpr int dimIntersection = 2;
76};
77
79template<class ct, int dw>
81{
82 using ctype = ct;
83 using Point = Dune::FieldVector<ctype, dw>;
84 using Polyhedron = std::vector<Point>;
85
86 static constexpr int dimWorld = dw;
87 static constexpr int dimIntersection = 3;
89};
90
92template<class Geometry1, class Geometry2>
94{
95 static constexpr int dimworld = Geometry1::coorddimension;
96 static constexpr int isDim = std::min( int(Geometry1::mydimension), int(Geometry2::mydimension) );
97 static_assert(dimworld == int(Geometry2::coorddimension),
98 "Geometries must have the same coordinate dimension!");
99 static_assert(int(Geometry1::mydimension) <= 3 && int(Geometry2::mydimension) <= 3,
100 "Geometries must have dimension 3 or less.");
101 using ctype = typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
102
103 using DefaultPolicies = std::tuple<PointPolicy<ctype, dimworld>,
107public:
108 using type = std::tuple_element_t<isDim, DefaultPolicies>;
109};
110
112template<class Geometry1, class Geometry2>
114
115} // end namespace IntersectionPolicy
116
117namespace Detail {
118
132template< class Geo1, class Geo2, class ctype,
133 class GetFacetCornerIndices, class ComputeNormalFunction >
134bool computeSegmentIntersection(const Geo1& geo1, const Geo2& geo2, ctype baseEps,
135 ctype& tfirst, ctype& tlast,
136 const GetFacetCornerIndices& getFacetCornerIndices,
137 const ComputeNormalFunction& computeNormal)
138{
139 static_assert(int(Geo2::mydimension) == 1, "Geometry2 must be a segment");
140 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
141 "Dimension of Geometry1 must be higher than that of Geometry2");
142
143 const auto a = geo2.corner(0);
144 const auto b = geo2.corner(1);
145 const auto d = b - a;
146
147 // The initial interval is the whole segment.
148 // Afterwards we start clipping the interval
149 // at the intersections with the facets of geo1
150 tfirst = 0.0;
151 tlast = 1.0;
152
153 const auto facets = getFacetCornerIndices(geo1);
154 for (const auto& f : facets)
155 {
156 const auto n = computeNormal(f);
157
158 const auto c0 = geo1.corner(f[0]);
159 const ctype denom = n*d;
160 const ctype dist = n*(a-c0);
161
162 // use first edge of the facet to scale eps
163 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
164 const ctype eps = baseEps*edge1.two_norm();
165
166 // if denominator is zero the segment in parallel to the edge.
167 // If the distance is positive there is no intersection
168 using std::abs;
169 if (abs(denom) < eps)
170 {
171 if (dist > eps)
172 return false;
173 }
174 else // not parallel: compute line-line intersection
175 {
176 const ctype t = -dist / denom;
177 // if entering half space cut tfirst if t is larger
178 using std::signbit;
179 if (signbit(denom))
180 {
181 if (t > tfirst)
182 tfirst = t;
183 }
184 // if exiting half space cut tlast if t is smaller
185 else
186 {
187 if (t < tlast)
188 tlast = t;
189 }
190 // there is no intersection of the interval is empty
191 // use unscaled epsilon since t is in local coordinates
192 if (tfirst > tlast - baseEps)
193 return false;
194 }
195 }
196
197 // If we made it until here an intersection exists.
198 return true;
199}
200
201} // end namespace Detail
202
209template
210<class Geometry1, class Geometry2,
211 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
212 int dimworld = Geometry1::coorddimension,
213 int dim1 = Geometry1::mydimension,
214 int dim2 = Geometry2::mydimension>
216{
217 static constexpr int dimWorld = Policy::dimWorld;
218
219public:
220 using ctype = typename Policy::ctype;
221 using Point = typename Policy::Point;
222 using Intersection = typename Policy::Intersection;
223
225 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
226 {
227 static_assert(dimworld == Geometry2::coorddimension, "Can only intersect geometries of same coordinate dimension");
228 DUNE_THROW(Dune::NotImplemented, "Geometry intersection detection with intersection computation not implemented for dimworld = "
229 << dimworld << ", dim1 = " << dim1 << ", dim2 = " << dim2);
230 }
231};
232
237template <class Geometry1, class Geometry2, class Policy>
238class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 1>
239{
240 enum { dimworld = 2 };
241 enum { dim1 = 1 };
242 enum { dim2 = 1 };
243
244 // base epsilon for floating point comparisons
245 static constexpr typename Policy::ctype eps_ = 1.5e-7;
246
247public:
248 using ctype = typename Policy::ctype;
249 using Point = typename Policy::Point;
250 using Intersection = typename Policy::Intersection;
251
258 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
259 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
260 {
261 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
262
263 const auto v1 = geo1.corner(1) - geo1.corner(0);
264 const auto v2 = geo2.corner(1) - geo2.corner(0);
265 const auto ac = geo2.corner(0) - geo1.corner(0);
266
267 auto n2 = Point({-1.0*v2[1], v2[0]});
268 n2 /= n2.two_norm();
269
270 // compute distance of first corner on geo2 to line1
271 const auto dist12 = n2*ac;
272
273 // first check parallel segments
274 using std::abs;
275 using std::sqrt;
276
277 const auto v1Norm2 = v1.two_norm2();
278 const auto eps = eps_*sqrt(v1Norm2);
279 const auto eps2 = eps_*v1Norm2;
280
281 const auto sp = n2*v1;
282 if (abs(sp) < eps)
283 {
284 if (abs(dist12) > eps)
285 return false;
286
287 // point intersection can only be one of the corners
288 if ( (v1*v2) < 0.0 )
289 {
290 if ( ac.two_norm2() < eps2 )
291 { intersection = geo2.corner(0); return true; }
292
293 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
294 { intersection = geo2.corner(1); return true; }
295 }
296 else
297 {
298 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
299 { intersection = geo2.corner(1); return true; }
300
301 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
302 { intersection = geo2.corner(0); return true; }
303 }
304
305 // no intersection
306 return false;
307 }
308
309 // intersection point on v1 in local coords
310 const auto t1 = dist12 / sp;
311
312 // check if the local coords are valid
313 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
314 return false;
315
316 // compute global coordinates
317 auto isPoint = geo1.global(t1);
318
319 // check if point is in bounding box of v2
320 const auto c = geo2.corner(0);
321 const auto d = geo2.corner(1);
322
323 using std::min; using std::max;
324 std::array<ctype, 4> bBox({ min(c[0], d[0]), min(c[1], d[1]), max(c[0], d[0]), max(c[1], d[1]) });
325
326 if ( intersectsPointBoundingBox(isPoint, bBox.data()) )
327 {
328 intersection = std::move(isPoint);
329 return true;
330 }
331
332 return false;
333 }
334
341 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
342 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
343 {
344 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
345
346 const auto& a = geo1.corner(0);
347 const auto& b = geo1.corner(1);
348 const auto ab = b - a;
349
350 const auto& p = geo2.corner(0);
351 const auto& q = geo2.corner(1);
352 const auto pq = q - p;
353 Dune::FieldVector<ctype, dimworld> n2{-pq[1], pq[0]};
354
355 using std::max;
356 const auto abNorm2 = ab.two_norm2();
357 const auto pqNorm2 = pq.two_norm2();
358 const auto eps2 = eps_*max(abNorm2, pqNorm2);
359
360 // non-parallel segments do not intersect in a segment.
361 using std::abs;
362 if (abs(n2*ab) > eps2)
363 return false;
364
365 // check if the segments lie on the same line
366 const auto ap = p - a;
367 if (abs(ap*n2) > eps2)
368 return false;
369
370 // compute scaled local coordinates of corner 1/2 of segment2 on segment1
371 auto t1 = ab*ap;
372 auto t2 = ab*(q - a);
373
374 using std::swap;
375 if (t1 > t2)
376 swap(t1, t2);
377
378 using std::clamp;
379 t1 = clamp(t1, 0.0, abNorm2);
380 t2 = clamp(t2, 0.0, abNorm2);
381
382 if (abs(t2-t1) < eps2)
383 return false;
384
385 intersection = Intersection({geo1.global(t1/abNorm2), geo1.global(t2/abNorm2)});
386 return true;
387 }
388};
389
394template <class Geometry1, class Geometry2, class Policy>
395class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 1>
396{
397 enum { dimworld = 2 };
398 enum { dim1 = 2 };
399 enum { dim2 = 1 };
400
401public:
402 using ctype = typename Policy::ctype;
403 using Point = typename Policy::Point;
404 using Intersection = typename Policy::Intersection;
405
406private:
407 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
408
409public:
417 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
418 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
419 {
420 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
421
422 ctype tfirst, tlast;
423 if (intersect_(geo1, geo2, tfirst, tlast))
424 {
425 // the intersection exists. Export the intersections geometry now:
426 // s(t) = a + t(b-a) in [tfirst, tlast]
427 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
428 return true;
429 }
430
431 return false;
432 }
433
441 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
442 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename P::Intersection& intersection)
443 {
444 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
445
446 ctype tfirst, tlast;
447 if (!intersect_(geo1, geo2, tfirst, tlast))
448 {
449 // if start & end point are same, the intersection is a point
450 using std::abs;
451 if ( tfirst > tlast - eps_ )
452 {
453 intersection = Intersection(geo2.global(tfirst));
454 return true;
455 }
456 }
457
458 return false;
459 }
460
461private:
465 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
466 {
467 // lambda to obtain the facet corners on geo1
468 auto getFacetCorners = [] (const Geometry1& geo1)
469 {
470 std::vector< std::array<int, 2> > facetCorners;
471 switch (geo1.corners())
472 {
473 case 4: // quadrilateral
474 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
475 break;
476 case 3: // triangle
477 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
478 break;
479 default:
480 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
481 << geo1.type() << " with "<< geo1.corners() << " corners.");
482 }
483
484 return facetCorners;
485 };
486
487 // lambda to obtain the normal on a facet
488 const auto center1 = geo1.center();
489 auto computeNormal = [&geo1, &center1] (const std::array<int, 2>& facetCorners)
490 {
491 const auto c0 = geo1.corner(facetCorners[0]);
492 const auto c1 = geo1.corner(facetCorners[1]);
493 const auto edge = c1 - c0;
494
495 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
496 n /= n.two_norm();
497
498 // orientation of the element is not known
499 // make sure normal points outwards of element
500 if ( n*(center1-c0) > 0.0 )
501 n *= -1.0;
502
503 return n;
504 };
505
506 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
507 }
508};
509
514template <class Geometry1, class Geometry2, class Policy>
515class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 2>
516: public GeometryIntersection<Geometry2, Geometry1, Policy, 2, 2, 1>
517{
519
520public:
528 template<class P = Policy>
529 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
530 {
531 return Base::intersection(geo2, geo1, intersection);
532 }
533};
534
539template <class Geometry1, class Geometry2, class Policy>
540class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 2>
541{
542 enum { dimworld = 2 };
543 enum { dim1 = 2 };
544 enum { dim2 = 2 };
545
546public:
547 using ctype = typename Policy::ctype;
548 using Point = typename Policy::Point;
549 using Intersection = typename Policy::Intersection;
550
551private:
552 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
553
554public:
567 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
568 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
569 {
570 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
571
572 // the candidate intersection points
573 std::vector<Point> points; points.reserve(6);
574
575 // add polygon1 corners that are inside polygon2
576 for (int i = 0; i < geo1.corners(); ++i)
577 if (intersectsPointGeometry(geo1.corner(i), geo2))
578 points.emplace_back(geo1.corner(i));
579
580 const auto numPoints1 = points.size();
581 if (numPoints1 != geo1.corners())
582 {
583 // add polygon2 corners that are inside polygon1
584 for (int i = 0; i < geo2.corners(); ++i)
585 if (intersectsPointGeometry(geo2.corner(i), geo1))
586 points.emplace_back(geo2.corner(i));
587
588 if (points.empty())
589 return false;
590
591 if (points.size() - numPoints1 != geo2.corners())
592 {
593 const auto refElement1 = referenceElement(geo1);
594 const auto refElement2 = referenceElement(geo2);
595
596 // add intersections of edges
597 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
599 for (int i = 0; i < refElement1.size(dim1-1); ++i)
600 {
601 const auto localEdgeGeom1 = refElement1.template geometry<dim1-1>(i);
602 const auto edge1 = SegGeometry( Dune::GeometryTypes::line,
603 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
604 geo1.global(localEdgeGeom1.corner(1))} ));
605
606 for (int j = 0; j < refElement2.size(dim2-1); ++j)
607 {
608 const auto localEdgeGeom2 = refElement2.template geometry<dim2-1>(j);
609 const auto edge2 = SegGeometry( Dune::GeometryTypes::line,
610 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
611 geo2.global(localEdgeGeom2.corner(1))} ));
612
614 typename EdgeTest::Intersection edgeIntersection;
615 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
616 points.emplace_back(edgeIntersection);
617 }
618 }
619 }
620 }
621
622 if (points.empty())
623 return false;
624
625 // remove duplicates
626 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
627 std::sort(points.begin(), points.end(), [&eps](const auto& a, const auto& b) -> bool
628 {
629 using std::abs;
630 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
631 });
632
633 auto removeIt = std::unique(points.begin(), points.end(), [&eps](const auto& a, const auto&b)
634 {
635 return (b-a).two_norm() < eps;
636 });
637
638 points.erase(removeIt, points.end());
639
640 // return false if we don't have at least three unique points
641 if (points.size() < 3)
642 return false;
643
644 // intersection polygon is convex hull of above points
645 intersection = grahamConvexHull<2>(points);
646 assert(!intersection.empty());
647 return true;
648 }
649
657 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
658 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
659 {
660 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
661 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for segment-like intersections");
662 }
663
671 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
672 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
673 {
674 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
675 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for touching points");
676 }
677};
678
683template <class Geometry1, class Geometry2, class Policy>
684class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 1>
685{
686 enum { dimworld = 3 };
687 enum { dim1 = 3 };
688 enum { dim2 = 1 };
689
690public:
691 using ctype = typename Policy::ctype;
692 using Point = typename Policy::Point;
693 using Intersection = typename Policy::Intersection;
694
695private:
696 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
697
698public:
710 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
711 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
712 {
713 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
714
715 ctype tfirst, tlast;
716 if (intersect_(geo1, geo2, tfirst, tlast))
717 {
718 // Intersection exists. Export the intersections geometry now:
719 // s(t) = a + t(b-a) in [tfirst, tlast]
720 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
721 return true;
722 }
723
724 return false;
725 }
726
738 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
739 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
740 {
741 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
742
743 ctype tfirst, tlast;
744 if (intersect_(geo1, geo2, tfirst, tlast))
745 {
746 // if start & end point are same, the intersection is a point
747 using std::abs;
748 if ( abs(tfirst - tlast) < eps_ )
749 {
750 intersection = Intersection(geo2.global(tfirst));
751 return true;
752 }
753 }
754
755 return false;
756 }
757
758private:
762 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
763 {
764 // lambda to obtain facet corners on geo1
765 auto getFacetCorners = [] (const Geometry1& geo1)
766 {
767 std::vector< std::vector<int> > facetCorners;
768 // sort facet corner so that normal n = (p1-p0)x(p2-p0) always points outwards
769 switch (geo1.corners())
770 {
771 // todo compute intersection geometries
772 case 8: // hexahedron
773 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
774 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
775 break;
776 case 4: // tetrahedron
777 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
778 break;
779 default:
780 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
781 << geo1.type() << ", "<< geo1.corners() << " corners.");
782 }
783
784 return facetCorners;
785 };
786
787 // lambda to obtain the normal on a facet
788 auto computeNormal = [&geo1] (const std::vector<int>& facetCorners)
789 {
790 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
791 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
792
793 auto n = crossProduct(v0, v1);
794 n /= n.two_norm();
795
796 return n;
797 };
798
799 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
800 }
801};
802
807template <class Geometry1, class Geometry2, class Policy>
808class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 3>
809: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 1>
810{
812public:
820 template<class P = Policy>
821 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
822 {
823 return Base::intersection(geo2, geo1, intersection);
824 }
825};
826
831template <class Geometry1, class Geometry2, class Policy>
832class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 2>
833{
834 enum { dimworld = 3 };
835 enum { dim1 = 3 };
836 enum { dim2 = 2 };
837
838public:
839 using ctype = typename Policy::ctype;
840 using Point = typename Policy::Point;
841 using Intersection = typename Policy::Intersection;
842
843private:
844 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
845
846public:
861 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
862 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
863 {
864 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
865
866 // the candidate intersection points
867 std::vector<Point> points; points.reserve(10);
868
869 // add 3d geometry corners that are inside the 2d geometry
870 for (int i = 0; i < geo1.corners(); ++i)
871 if (intersectsPointGeometry(geo1.corner(i), geo2))
872 points.emplace_back(geo1.corner(i));
873
874 // add 2d geometry corners that are inside the 3d geometry
875 for (int i = 0; i < geo2.corners(); ++i)
876 if (intersectsPointGeometry(geo2.corner(i), geo1))
877 points.emplace_back(geo2.corner(i));
878
879 // get some geometry types
880 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
881 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
882
883 const auto refElement1 = referenceElement(geo1);
884 const auto refElement2 = referenceElement(geo2);
885
886 // intersection policy for point-like intersections (used later)
888
889 // add intersection points of all polyhedron edges (codim dim-1) with the polygon
890 for (int i = 0; i < refElement1.size(dim1-1); ++i)
891 {
892 const auto localEdgeGeom = refElement1.template geometry<dim1-1>(i);
893 const auto p = geo1.global(localEdgeGeom.corner(0));
894 const auto q = geo1.global(localEdgeGeom.corner(1));
895 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
896
898 typename PolySegTest::Intersection polySegIntersection;
899 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
900 points.emplace_back(polySegIntersection);
901 }
902
903 // add intersection points of all polygon faces (codim 1) with the polyhedron faces
904 for (int i = 0; i < refElement1.size(1); ++i)
905 {
906 const auto faceGeo = [&]()
907 {
908 const auto localFaceGeo = refElement1.template geometry<1>(i);
909 if (localFaceGeo.corners() == 4)
910 {
911 const auto a = geo1.global(localFaceGeo.corner(0));
912 const auto b = geo1.global(localFaceGeo.corner(1));
913 const auto c = geo1.global(localFaceGeo.corner(2));
914 const auto d = geo1.global(localFaceGeo.corner(3));
915
916 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
917 }
918 else
919 {
920 const auto a = geo1.global(localFaceGeo.corner(0));
921 const auto b = geo1.global(localFaceGeo.corner(1));
922 const auto c = geo1.global(localFaceGeo.corner(2));
923
924 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
925 }
926 }();
927
928 for (int j = 0; j < refElement2.size(1); ++j)
929 {
930 const auto localEdgeGeom = refElement2.template geometry<1>(j);
931 const auto p = geo2.global(localEdgeGeom.corner(0));
932 const auto q = geo2.global(localEdgeGeom.corner(1));
933
934 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
935
937 typename PolySegTest::Intersection polySegIntersection;
938 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
939 points.emplace_back(polySegIntersection);
940 }
941 }
942
943 // return if no intersection points were found
944 if (points.empty()) return false;
945
946 // remove duplicates
947 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
948 std::sort(points.begin(), points.end(), [&eps](const auto& a, const auto& b) -> bool
949 {
950 using std::abs;
951 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : (abs(a[1]-b[1]) > eps ? a[1] < b[1] : (a[2] < b[2])));
952 });
953
954 auto removeIt = std::unique(points.begin(), points.end(), [&eps](const auto& a, const auto&b)
955 {
956 return (b-a).two_norm() < eps;
957 });
958
959 points.erase(removeIt, points.end());
960
961 // return false if we don't have more than three unique points
962 if (points.size() < 3) return false;
963
964 // intersection polygon is convex hull of above points
965 intersection = grahamConvexHull<2>(points);
966 assert(!intersection.empty());
967
968 return true;
969 }
970
985 template<class P = Policy, std::enable_if_t<P::dimIntersection != 2, int> = 0>
986 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
987 {
988 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
989 DUNE_THROW(Dune::NotImplemented, "Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
990 }
991};
992
997template <class Geometry1, class Geometry2, class Policy>
998class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 3>
999: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 2>
1000{
1002public:
1010 template<class P = Policy>
1011 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
1012 {
1013 return Base::intersection(geo2, geo1, intersection);
1014 }
1015};
1016
1021template <class Geometry1, class Geometry2, class Policy>
1022class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 1>
1023{
1024 enum { dimworld = 3 };
1025 enum { dim1 = 2 };
1026 enum { dim2 = 1 };
1027
1028public:
1029 using ctype = typename Policy::ctype;
1030 using Point = typename Policy::Point;
1031 using Intersection = typename Policy::Intersection;
1032
1033private:
1034 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1035
1036public:
1046 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1047 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1048 {
1049 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1050
1051 ctype tfirst, tlast;
1052 if (intersect_(geo1, geo2, tfirst, tlast))
1053 {
1054 // the intersection exists. Export the intersections geometry now:
1055 // s(t) = a + t(b-a) in [tfirst, tlast]
1056 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
1057 return true;
1058 }
1059
1060 return false;
1061 }
1062
1072 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1073 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& is)
1074 {
1075 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1076
1077 const auto p = geo2.corner(0);
1078 const auto q = geo2.corner(1);
1079
1080 const auto a = geo1.corner(0);
1081 const auto b = geo1.corner(1);
1082 const auto c = geo1.corner(2);
1083
1084 if (geo1.corners() == 3)
1085 return intersect_<Policy>(a, b, c, p, q, is);
1086
1087 else if (geo1.corners() == 4)
1088 {
1089 Intersection is1, is2;
1090 bool hasSegment1, hasSegment2;
1091
1092 const auto d = geo1.corner(3);
1093 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1094 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1095
1096 if (hasSegment1 || hasSegment2)
1097 return false;
1098
1099 if (intersects1) { is = is1; return true; }
1100 if (intersects2) { is = is2; return true; }
1101
1102 return false;
1103 }
1104
1105 else
1106 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1107 << geo1.type() << ", "<< geo1.corners() << " corners.");
1108 }
1109
1110private:
1114 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1115 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
1116 {
1117 // lambda to obtain the facet corners on geo1
1118 auto getFacetCorners = [] (const Geometry1& geo1)
1119 {
1120 std::vector< std::array<int, 2> > facetCorners;
1121 switch (geo1.corners())
1122 {
1123 case 4: // quadrilateral
1124 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1125 break;
1126 case 3: // triangle
1127 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1128 break;
1129 default:
1130 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1131 << geo1.type() << " with "<< geo1.corners() << " corners.");
1132 }
1133
1134 return facetCorners;
1135 };
1136
1137 const auto center1 = geo1.center();
1138 const auto normal1 = crossProduct(geo1.corner(1) - geo1.corner(0),
1139 geo1.corner(2) - geo1.corner(0));
1140
1141 // lambda to obtain the normal on a facet
1142 auto computeNormal = [&center1, &normal1, &geo1] (const std::array<int, 2>& facetCorners)
1143 {
1144 const auto c0 = geo1.corner(facetCorners[0]);
1145 const auto c1 = geo1.corner(facetCorners[1]);
1146 const auto edge = c1 - c0;
1147
1148 Dune::FieldVector<ctype, dimworld> n = crossProduct(edge, normal1);
1149 n /= n.two_norm();
1150
1151 // orientation of the element is not known
1152 // make sure normal points outwards of element
1153 if ( n*(center1-c0) > 0.0 )
1154 n *= -1.0;
1155
1156 return n;
1157 };
1158
1159 return Detail::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
1160 }
1161
1165 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1166 static bool intersect_(const Point& a, const Point& b, const Point& c,
1167 const Point& p, const Point& q,
1168 Intersection& is)
1169 {
1170 bool hasSegment;
1171 return intersect_(a, b, c, p, q, is, hasSegment);
1172 }
1173
1178 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1179 static bool intersect_(const Point& a, const Point& b, const Point& c,
1180 const Point& p, const Point& q,
1181 Intersection& is, bool& hasSegmentIs)
1182 {
1183 hasSegmentIs = false;
1184
1185 const auto ab = b - a;
1186 const auto ac = c - a;
1187 const auto qp = p - q;
1188
1189 // compute the triangle normal that defines the triangle plane
1190 const auto normal = crossProduct(ab, ac);
1191
1192 // compute the denominator
1193 // if denom is 0 the segment is parallel and we can return
1194 const auto denom = normal*qp;
1195 const auto abnorm2 = ab.two_norm2();
1196 const auto eps = eps_*abnorm2*qp.two_norm();
1197 using std::abs;
1198 if (abs(denom) < eps)
1199 {
1200 const auto pa = a - p;
1201 const auto denom2 = normal*pa;
1202 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1203 return false;
1204
1205 // if we get here, we are in-plane. Check if a
1206 // segment-like intersection with geometry 1 exists.
1207 using SegmentPolicy = typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1208 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1209 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1210 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1211 using SegmentIntersectionType = typename SegmentIntersectionAlgorithm::Intersection;
1212 SegmentIntersectionType segmentIs;
1213
1214 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1215 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1216 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1217 {
1218 hasSegmentIs = true;
1219 return false;
1220 }
1221
1222 // We are in-plane. A point-like
1223 // intersection can only be on the edges
1224 for (const auto& ip : {p, q})
1225 {
1226 if ( intersectsPointSimplex(ip, a, b)
1227 || intersectsPointSimplex(ip, b, c)
1228 || intersectsPointSimplex(ip, c, a) )
1229 {
1230 is = ip;
1231 return true;
1232 }
1233 }
1234
1235 return false;
1236 }
1237
1238 // compute intersection t value of pq with plane of triangle.
1239 // a segment intersects if and only if 0 <= t <= 1.
1240 const auto ap = p - a;
1241 const auto t = (ap*normal)/denom;
1242 if (t < 0.0 - eps_) return false;
1243 if (t > 1.0 + eps_) return false;
1244
1245 // compute the barycentric coordinates and check if the intersection point
1246 // is inside the bounds of the triangle
1247 const auto e = crossProduct(qp, ap);
1248 const auto v = (ac*e)/denom;
1249 if (v < -eps_ || v > 1.0 + eps_) return false;
1250 const auto w = -(ab*e)/denom;
1251 if (w < -eps_ || v + w > 1.0 + eps_) return false;
1252
1253 // Now we are sure there is an intersection points
1254 // Perform delayed division compute the last barycentric coordinate component
1255 const auto u = 1.0 - v - w;
1256
1257 Point ip(0.0);
1258 ip.axpy(u, a);
1259 ip.axpy(v, b);
1260 ip.axpy(w, c);
1261 is = ip;
1262 return true;
1263 }
1264};
1265
1270template <class Geometry1, class Geometry2, class Policy>
1271class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 2>
1272: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 2, 1>
1273{
1275public:
1283 template<class P = Policy>
1284 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
1285 {
1286 return Base::intersection(geo2, geo1, intersection);
1287 }
1288};
1289
1294template <class Geometry1, class Geometry2, class Policy>
1295class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 1>
1296{
1297 enum { dimworld = 3 };
1298 enum { dim1 = 1 };
1299 enum { dim2 = 1 };
1300
1301public:
1302 using ctype = typename Policy::ctype;
1303 using Point = typename Policy::Point;
1304 using Intersection = typename Policy::Intersection;
1305
1306private:
1307 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1308
1309public:
1317 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1318 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1319 {
1320 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1321 DUNE_THROW(Dune::NotImplemented, "segment-segment intersection detection for point-like intersections");
1322 }
1323
1331 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1332 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1333 {
1334 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1335
1336 const auto& a = geo1.corner(0);
1337 const auto& b = geo1.corner(1);
1338 const auto ab = b-a;
1339
1340 const auto& p = geo2.corner(0);
1341 const auto& q = geo2.corner(1);
1342 const auto pq = q-p;
1343
1344 const auto abNorm2 = ab.two_norm2();
1345 const auto pqNorm2 = pq.two_norm2();
1346
1347 using std::max;
1348 const auto eps2 = eps_*max(abNorm2, pqNorm2);
1349
1350 // if the segment are not parallel there is no segment intersection
1351 using std::abs;
1352 if (crossProduct(ab, pq).two_norm2() > eps2*eps2)
1353 return false;
1354
1355 const auto ap = (p-a);
1356 const auto aq = (q-a);
1357
1358 // points have to be colinear
1359 if (crossProduct(ap, aq).two_norm2() > eps2*eps2)
1360 return false;
1361
1362 // scaled local coordinates
1363 // we save the division for the normalization for now
1364 // and do it in the very end if we find an intersection
1365 auto tp = ap*ab;
1366 auto tq = aq*ab;
1367
1368 // make sure they are sorted
1369 using std::swap;
1370 if (tp > tq)
1371 swap(tp, tq);
1372
1373 using std::clamp;
1374 tp = clamp(tp, 0.0, abNorm2);
1375 tq = clamp(tq, 0.0, abNorm2);
1376
1377 if (abs(tp-tq) < eps2)
1378 return false;
1379
1380 intersection = Intersection({geo1.global(tp/abNorm2), geo1.global(tq/abNorm2)});
1381 return true;
1382 }
1383};
1384
1385} // end namespace Dumux
1386
1387#endif
Define some often used mathematical functions.
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: geometry/intersectspointgeometry.hh:38
Vector normal(const Vector &v)
Create a vector normal to the given one (v is expected to be non-zero)
Definition: geometry/normal.hh:36
bool intersectsPointSimplex(const Dune::FieldVector< ctype, dimworld > &point, const Dune::FieldVector< ctype, dimworld > &p0, const Dune::FieldVector< ctype, dimworld > &p1, const Dune::FieldVector< ctype, dimworld > &p2, const Dune::FieldVector< ctype, dimworld > &p3)
Find out whether a point is inside a tetrahedron (p0, p1, p2, p3) (dimworld is 3)
Definition: geometry/intersectspointsimplex.hh:36
bool computeSegmentIntersection(const Geo1 &geo1, const Geo2 &geo2, ctype baseEps, ctype &tfirst, ctype &tlast, const GetFacetCornerIndices &getFacetCornerIndices, const ComputeNormalFunction &computeNormal)
Algorithm to find segment-like intersections of a polgon/polyhedron with a segment....
Definition: geometry/geometryintersection.hh:134
Dune::FieldVector< Scalar, 3 > crossProduct(const Dune::FieldVector< Scalar, 3 > &vec1, const Dune::FieldVector< Scalar, 3 > &vec2)
Cross product of two vectors in three-dimensional Euclidean space.
Definition: math.hh:640
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: geometry/boundingboxtree.hh:304
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometry/geometryintersection.hh:113
constexpr Line line
Definition: couplingmanager1d3d_line.hh:52
Policy structure for point-like intersections.
Definition: geometry/geometryintersection.hh:43
static constexpr int dimIntersection
Definition: geometry/geometryintersection.hh:48
Point Intersection
Definition: geometry/geometryintersection.hh:49
ct ctype
Definition: geometry/geometryintersection.hh:44
Dune::FieldVector< ctype, dw > Point
Definition: geometry/geometryintersection.hh:45
static constexpr int dimWorld
Definition: geometry/geometryintersection.hh:47
Policy structure for segment-like intersections.
Definition: geometry/geometryintersection.hh:55
Segment Intersection
Definition: geometry/geometryintersection.hh:62
static constexpr int dimIntersection
Definition: geometry/geometryintersection.hh:61
Dune::FieldVector< ctype, dw > Point
Definition: geometry/geometryintersection.hh:57
static constexpr int dimWorld
Definition: geometry/geometryintersection.hh:60
std::array< Point, 2 > Segment
Definition: geometry/geometryintersection.hh:58
ct ctype
Definition: geometry/geometryintersection.hh:56
Policy structure for polygon-like intersections.
Definition: geometry/geometryintersection.hh:68
static constexpr int dimWorld
Definition: geometry/geometryintersection.hh:73
static constexpr int dimIntersection
Definition: geometry/geometryintersection.hh:74
ct ctype
Definition: geometry/geometryintersection.hh:69
Polygon Intersection
Definition: geometry/geometryintersection.hh:75
Dune::FieldVector< ctype, dw > Point
Definition: geometry/geometryintersection.hh:70
std::vector< Point > Polygon
Definition: geometry/geometryintersection.hh:71
Policy structure for polyhedron-like intersections.
Definition: geometry/geometryintersection.hh:81
std::vector< Point > Polyhedron
Definition: geometry/geometryintersection.hh:84
Dune::FieldVector< ctype, dw > Point
Definition: geometry/geometryintersection.hh:83
Polyhedron Intersection
Definition: geometry/geometryintersection.hh:88
static constexpr int dimIntersection
Definition: geometry/geometryintersection.hh:87
static constexpr int dimWorld
Definition: geometry/geometryintersection.hh:86
ct ctype
Definition: geometry/geometryintersection.hh:82
default policy chooser class
Definition: geometry/geometryintersection.hh:94
std::tuple_element_t< isDim, DefaultPolicies > type
Definition: geometry/geometryintersection.hh:108
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometry/geometryintersection.hh:216
typename Policy::ctype ctype
Definition: geometry/geometryintersection.hh:220
typename Policy::Intersection Intersection
Definition: geometry/geometryintersection.hh:222
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Determine if the two geometries intersect and compute the intersection geometry.
Definition: geometry/geometryintersection.hh:225
typename Policy::Point Point
Definition: geometry/geometryintersection.hh:221
typename Policy::Intersection Intersection
Definition: geometry/geometryintersection.hh:250
typename Policy::ctype ctype
Definition: geometry/geometryintersection.hh:248
typename Policy::Point Point
Definition: geometry/geometryintersection.hh:249
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometry/geometryintersection.hh:259
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometry/geometryintersection.hh:418
typename Policy::ctype ctype
Definition: geometry/geometryintersection.hh:402
typename Policy::Intersection Intersection
Definition: geometry/geometryintersection.hh:404
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometry/geometryintersection.hh:442
typename Policy::Point Point
Definition: geometry/geometryintersection.hh:403
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometry/geometryintersection.hh:529
typename Policy::Intersection Intersection
Definition: geometry/geometryintersection.hh:549
typename Policy::ctype ctype
Definition: geometry/geometryintersection.hh:547
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometry/geometryintersection.hh:568
typename Policy::Point Point
Definition: geometry/geometryintersection.hh:548
typename Policy::Intersection Intersection
Definition: geometry/geometryintersection.hh:693
typename Policy::Point Point
Definition: geometry/geometryintersection.hh:692
typename Policy::ctype ctype
Definition: geometry/geometryintersection.hh:691
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometry/geometryintersection.hh:711
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometry/geometryintersection.hh:821
typename Policy::Intersection Intersection
Definition: geometry/geometryintersection.hh:841
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometry/geometryintersection.hh:862
typename Policy::Point Point
Definition: geometry/geometryintersection.hh:840
typename Policy::ctype ctype
Definition: geometry/geometryintersection.hh:839
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometry/geometryintersection.hh:1011
typename Policy::Intersection Intersection
Definition: geometry/geometryintersection.hh:1031
typename Policy::Point Point
Definition: geometry/geometryintersection.hh:1030
typename Policy::ctype ctype
Definition: geometry/geometryintersection.hh:1029
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometry/geometryintersection.hh:1047
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometry/geometryintersection.hh:1073
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometry/geometryintersection.hh:1284
typename Policy::Intersection Intersection
Definition: geometry/geometryintersection.hh:1304
typename Policy::Point Point
Definition: geometry/geometryintersection.hh:1303
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometry/geometryintersection.hh:1318
typename Policy::ctype ctype
Definition: geometry/geometryintersection.hh:1302
An axis-aligned bounding box volume hierarchy for dune grids.
A function to compute the convex hull of a point cloud and a function to triangulate the polygon area...
Detect if a point intersects a geometry.