3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
geometryintersection.hh
Go to the documentation of this file.
1/*****************************************************************************
2 * See the file COPYING for full copying permissions. *
3 * *
4 * This program is free software: you can redistribute it and/or modify *
5 * it under the terms of the GNU General Public License as published by *
6 * the Free Software Foundation, either version 3 of the License, or *
7 * (at your option) any later version. *
8 * *
9 * This program is distributed in the hope that it will be useful, *
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
12 * GNU General Public License for more details. *
13 * *
14 * You should have received a copy of the GNU General Public License *
15 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
16 *****************************************************************************/
23#ifndef DUMUX_GEOMETRY_INTERSECTION_HH
24#define DUMUX_GEOMETRY_INTERSECTION_HH
25
26#include <tuple>
27
28#include <dune/common/exceptions.hh>
29#include <dune/common/promotiontraits.hh>
30#include <dune/geometry/referenceelements.hh>
31#include <dune/geometry/multilineargeometry.hh>
32
33#include <dumux/common/math.hh>
37
38namespace Dumux {
39namespace IntersectionPolicy {
40
42template<class ct, int dw>
44{
45 using ctype = ct;
46 using Point = Dune::FieldVector<ctype, dw>;
47
48 static constexpr int dimWorld = dw;
49 static constexpr int dimIntersection = 0;
51};
52
54template<class ct, int dw>
56{
57 using ctype = ct;
58 using Point = Dune::FieldVector<ctype, dw>;
59 using Segment = std::array<Point, 2>;
60
61 static constexpr int dimWorld = dw;
62 static constexpr int dimIntersection = 1;
64};
65
67template<class ct, int dw>
69{
70 using ctype = ct;
71 using Point = Dune::FieldVector<ctype, dw>;
72 using Polygon = std::vector<Point>;
73
74 static constexpr int dimWorld = dw;
75 static constexpr int dimIntersection = 2;
77};
78
80template<class ct, int dw>
82{
83 using ctype = ct;
84 using Point = Dune::FieldVector<ctype, dw>;
85 using Polyhedron = std::vector<Point>;
86
87 static constexpr int dimWorld = dw;
88 static constexpr int dimIntersection = 3;
90};
91
93template<class Geometry1, class Geometry2>
95{
96 static constexpr int dimworld = Geometry1::coorddimension;
97 static constexpr int isDim = std::min( int(Geometry1::mydimension), int(Geometry2::mydimension) );
98 static_assert(dimworld == int(Geometry2::coorddimension),
99 "Geometries must have the same coordinate dimension!");
100 static_assert(int(Geometry1::mydimension) <= 3 && int(Geometry2::mydimension) <= 3,
101 "Geometries must have dimension 3 or less.");
102 using ctype = typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
103
104 using DefaultPolicies = std::tuple<PointPolicy<ctype, dimworld>,
108public:
109 using type = std::tuple_element_t<isDim, DefaultPolicies>;
110};
111
113template<class Geometry1, class Geometry2>
115
116} // end namespace IntersectionPolicy
117
118namespace Impl {
119
133template< class Geo1, class Geo2, class ctype,
134 class GetFacetCornerIndices, class ComputeNormalFunction >
135bool computeSegmentIntersection(const Geo1& geo1, const Geo2& geo2, ctype baseEps,
136 ctype& tfirst, ctype& tlast,
137 const GetFacetCornerIndices& getFacetCornerIndices,
138 const ComputeNormalFunction& computeNormal)
139{
140 static_assert(int(Geo2::mydimension) == 1, "Geometry2 must be a segment");
141 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
142 "Dimension of Geometry1 must be higher than that of Geometry2");
143
144 const auto a = geo2.corner(0);
145 const auto b = geo2.corner(1);
146 const auto d = b - a;
147
148 // The initial interval is the whole segment.
149 // Afterwards we start clipping the interval
150 // at the intersections with the facets of geo1
151 tfirst = 0.0;
152 tlast = 1.0;
153
154 const auto facets = getFacetCornerIndices(geo1);
155 for (const auto& f : facets)
156 {
157 const auto n = computeNormal(f);
158
159 const auto c0 = geo1.corner(f[0]);
160 const ctype denom = n*d;
161 const ctype dist = n*(a-c0);
162
163 // use first edge of the facet to scale eps
164 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
165 const ctype eps = baseEps*edge1.two_norm();
166
167 // if denominator is zero the segment in parallel to the edge.
168 // If the distance is positive there is no intersection
169 using std::abs;
170 if (abs(denom) < eps)
171 {
172 if (dist > eps)
173 return false;
174 }
175 else // not parallel: compute line-line intersection
176 {
177 const ctype t = -dist / denom;
178 // if entering half space cut tfirst if t is larger
179 using std::signbit;
180 if (signbit(denom))
181 {
182 if (t > tfirst)
183 tfirst = t;
184 }
185 // if exiting half space cut tlast if t is smaller
186 else
187 {
188 if (t < tlast)
189 tlast = t;
190 }
191 // there is no intersection of the interval is empty
192 // use unscaled epsilon since t is in local coordinates
193 if (tfirst > tlast - baseEps)
194 return false;
195 }
196 }
197
198 // If we made it until here an intersection exists.
199 return true;
200}
201
202} // end namespace Impl
203
210template
211<class Geometry1, class Geometry2,
212 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
213 int dimworld = Geometry1::coorddimension,
214 int dim1 = Geometry1::mydimension,
215 int dim2 = Geometry2::mydimension>
217{
218 static constexpr int dimWorld = Policy::dimWorld;
219
220public:
221 using ctype = typename Policy::ctype;
222 using Point = typename Policy::Point;
223 using Intersection = typename Policy::Intersection;
224
226 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
227 {
228 static_assert(dimworld == Geometry2::coorddimension, "Can only intersect geometries of same coordinate dimension");
229 DUNE_THROW(Dune::NotImplemented, "Geometry intersection detection with intersection computation not implemented for dimworld = "
230 << dimworld << ", dim1 = " << dim1 << ", dim2 = " << dim2);
231 }
232};
233
238template <class Geometry1, class Geometry2, class Policy>
239class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 1>
240{
241 enum { dimworld = 2 };
242 enum { dim1 = 1 };
243 enum { dim2 = 1 };
244
245 // base epsilon for floating point comparisons
246 static constexpr typename Policy::ctype eps_ = 1.5e-7;
247
248public:
249 using ctype = typename Policy::ctype;
250 using Point = typename Policy::Point;
251 using Intersection = typename Policy::Intersection;
252
259 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
260 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
261 {
262 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
263
264 const auto v1 = geo1.corner(1) - geo1.corner(0);
265 const auto v2 = geo2.corner(1) - geo2.corner(0);
266 const auto ac = geo2.corner(0) - geo1.corner(0);
267
268 auto n2 = Point({-1.0*v2[1], v2[0]});
269 n2 /= n2.two_norm();
270
271 // compute distance of first corner on geo2 to line1
272 const auto dist12 = n2*ac;
273
274 // first check parallel segments
275 using std::abs;
276 using std::sqrt;
277
278 const auto v1Norm2 = v1.two_norm2();
279 const auto eps = eps_*sqrt(v1Norm2);
280 const auto eps2 = eps_*v1Norm2;
281
282 const auto sp = n2*v1;
283 if (abs(sp) < eps)
284 {
285 if (abs(dist12) > eps)
286 return false;
287
288 // point intersection can only be one of the corners
289 if ( (v1*v2) < 0.0 )
290 {
291 if ( ac.two_norm2() < eps2 )
292 { intersection = geo2.corner(0); return true; }
293
294 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
295 { intersection = geo2.corner(1); return true; }
296 }
297 else
298 {
299 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
300 { intersection = geo2.corner(1); return true; }
301
302 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
303 { intersection = geo2.corner(0); return true; }
304 }
305
306 // no intersection
307 return false;
308 }
309
310 // intersection point on v1 in local coords
311 const auto t1 = dist12 / sp;
312
313 // check if the local coords are valid
314 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
315 return false;
316
317 // compute global coordinates
318 auto isPoint = geo1.global(t1);
319
320 // check if point is in bounding box of v2
321 const auto c = geo2.corner(0);
322 const auto d = geo2.corner(1);
323
324 using std::min; using std::max;
325 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]) });
326
327 if ( intersectsPointBoundingBox(isPoint, bBox.data()) )
328 {
329 intersection = std::move(isPoint);
330 return true;
331 }
332
333 return false;
334 }
335
343 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
344 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
345 {
346 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
347 DUNE_THROW(Dune::NotImplemented, "segment-segment intersection detection for segment-like intersections");
348 }
349};
350
355template <class Geometry1, class Geometry2, class Policy>
356class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 1>
357{
358 enum { dimworld = 2 };
359 enum { dim1 = 2 };
360 enum { dim2 = 1 };
361
362public:
363 using ctype = typename Policy::ctype;
364 using Point = typename Policy::Point;
365 using Intersection = typename Policy::Intersection;
366
367private:
368 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
369 using ReferenceElements = typename Dune::ReferenceElements<ctype, dim1>;
370
371public:
379 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
380 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
381 {
382 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
383
384 ctype tfirst, tlast;
385 if (intersect_(geo1, geo2, tfirst, tlast))
386 {
387 // the intersection exists. Export the intersections geometry now:
388 // s(t) = a + t(b-a) in [tfirst, tlast]
389 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
390 return true;
391 }
392
393 return false;
394 }
395
403 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
404 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename P::Intersection& intersection)
405 {
406 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
407
408 ctype tfirst, tlast;
409 if (!intersect_(geo1, geo2, tfirst, tlast))
410 {
411 // if start & end point are same, the intersection is a point
412 using std::abs;
413 if ( tfirst > tlast - eps_ )
414 {
415 intersection = Intersection(geo2.global(tfirst));
416 return true;
417 }
418 }
419
420 return false;
421 }
422
423private:
427 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
428 {
429 // lambda to obtain the facet corners on geo1
430 auto getFacetCorners = [] (const Geometry1& geo1)
431 {
432 std::vector< std::array<int, 2> > facetCorners;
433 switch (geo1.corners())
434 {
435 case 4: // quadrilateral
436 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
437 break;
438 case 3: // triangle
439 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
440 break;
441 default:
442 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
443 << geo1.type() << " with "<< geo1.corners() << " corners.");
444 }
445
446 return facetCorners;
447 };
448
449 // lambda to obtain the normal on a facet
450 const auto center1 = geo1.center();
451 auto computeNormal = [&geo1, &center1] (const std::array<int, 2>& facetCorners)
452 {
453 const auto c0 = geo1.corner(facetCorners[0]);
454 const auto c1 = geo1.corner(facetCorners[1]);
455 const auto edge = c1 - c0;
456
457 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
458 n /= n.two_norm();
459
460 // orientation of the element is not known
461 // make sure normal points outwards of element
462 if ( n*(center1-c0) > 0.0 )
463 n *= -1.0;
464
465 return n;
466 };
467
468 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
469 }
470};
471
476template <class Geometry1, class Geometry2, class Policy>
477class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 2>
478: public GeometryIntersection<Geometry2, Geometry1, Policy, 2, 2, 1>
479{
481
482public:
490 template<class P = Policy>
491 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
492 {
493 return Base::intersection(geo2, geo1, intersection);
494 }
495};
496
501template <class Geometry1, class Geometry2, class Policy>
502class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 2>
503{
504 enum { dimworld = 2 };
505 enum { dim1 = 2 };
506 enum { dim2 = 2 };
507
508public:
509 using ctype = typename Policy::ctype;
510 using Point = typename Policy::Point;
511 using Intersection = typename Policy::Intersection;
512
513private:
514 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
515 using ReferenceElementsGeo1 = typename Dune::ReferenceElements<ctype, dim1>;
516 using ReferenceElementsGeo2 = typename Dune::ReferenceElements<ctype, dim2>;
517
518public:
531 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
532 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
533 {
534 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
535
536 // the candidate intersection points
537 std::vector<Point> points; points.reserve(6);
538
539 // add polygon1 corners that are inside polygon2
540 for (int i = 0; i < geo1.corners(); ++i)
541 if (intersectsPointGeometry(geo1.corner(i), geo2))
542 points.emplace_back(geo1.corner(i));
543
544 const auto numPoints1 = points.size();
545 if (numPoints1 != geo1.corners())
546 {
547 // add polygon2 corners that are inside polygon1
548 for (int i = 0; i < geo2.corners(); ++i)
549 if (intersectsPointGeometry(geo2.corner(i), geo1))
550 points.emplace_back(geo2.corner(i));
551
552 if (points.empty())
553 return false;
554
555 if (points.size() - numPoints1 != geo2.corners())
556 {
557 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
558 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
559
560 // add intersections of edges
561 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
563 for (int i = 0; i < referenceElement1.size(dim1-1); ++i)
564 {
565 const auto localEdgeGeom1 = referenceElement1.template geometry<dim1-1>(i);
566 const auto edge1 = SegGeometry( Dune::GeometryTypes::line,
567 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
568 geo1.global(localEdgeGeom1.corner(1))} ));
569
570 for (int j = 0; j < referenceElement2.size(dim2-1); ++j)
571 {
572 const auto localEdgeGeom2 = referenceElement2.template geometry<dim2-1>(j);
573 const auto edge2 = SegGeometry( Dune::GeometryTypes::line,
574 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
575 geo2.global(localEdgeGeom2.corner(1))} ));
576
578 typename EdgeTest::Intersection edgeIntersection;
579 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
580 points.emplace_back(edgeIntersection);
581 }
582 }
583 }
584 }
585
586 if (points.empty())
587 return false;
588
589 // remove duplicates
590 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
591 std::sort(points.begin(), points.end(), [&eps](const auto& a, const auto& b) -> bool
592 {
593 using std::abs;
594 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
595 });
596
597 auto removeIt = std::unique(points.begin(), points.end(), [&eps](const auto& a, const auto&b)
598 {
599 return (b-a).two_norm() < eps;
600 });
601
602 points.erase(removeIt, points.end());
603
604 // return false if we don't have at least three unique points
605 if (points.size() < 3)
606 return false;
607
608 // intersection polygon is convex hull of above points
609 intersection = grahamConvexHull<2>(points);
610 assert(!intersection.empty());
611 return true;
612 }
613
621 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
622 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
623 {
624 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
625 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for segment-like intersections");
626 }
627
635 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
636 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
637 {
638 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
639 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for touching points");
640 }
641};
642
647template <class Geometry1, class Geometry2, class Policy>
648class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 1>
649{
650 enum { dimworld = 3 };
651 enum { dim1 = 3 };
652 enum { dim2 = 1 };
653
654public:
655 using ctype = typename Policy::ctype;
656 using Point = typename Policy::Point;
657 using Intersection = typename Policy::Intersection;
658
659private:
660 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
661 using ReferenceElements = typename Dune::ReferenceElements<ctype, dim1>;
662
663public:
675 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
676 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
677 {
678 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
679
680 ctype tfirst, tlast;
681 if (intersect_(geo1, geo2, tfirst, tlast))
682 {
683 // Intersection exists. Export the intersections geometry now:
684 // s(t) = a + t(b-a) in [tfirst, tlast]
685 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
686 return true;
687 }
688
689 return false;
690 }
691
703 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
704 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
705 {
706 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
707
708 ctype tfirst, tlast;
709 if (intersect_(geo1, geo2, tfirst, tlast))
710 {
711 // if start & end point are same, the intersection is a point
712 using std::abs;
713 if ( abs(tfirst - tlast) < eps_ )
714 {
715 intersection = Intersection(geo2.global(tfirst));
716 return true;
717 }
718 }
719
720 return false;
721 }
722
723private:
727 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
728 {
729 // lambda to obtain facet corners on geo1
730 auto getFacetCorners = [] (const Geometry1& geo1)
731 {
732 std::vector< std::vector<int> > facetCorners;
733 // sort facet corner so that normal n = (p1-p0)x(p2-p0) always points outwards
734 switch (geo1.corners())
735 {
736 // todo compute intersection geometries
737 case 8: // hexahedron
738 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
739 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
740 break;
741 case 4: // tetrahedron
742 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
743 break;
744 default:
745 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
746 << geo1.type() << ", "<< geo1.corners() << " corners.");
747 }
748
749 return facetCorners;
750 };
751
752 // lambda to obtain the normal on a facet
753 auto computeNormal = [&geo1] (const std::vector<int>& facetCorners)
754 {
755 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
756 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
757
758 auto n = crossProduct(v0, v1);
759 n /= n.two_norm();
760
761 return n;
762 };
763
764 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
765 }
766};
767
772template <class Geometry1, class Geometry2, class Policy>
773class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 3>
774: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 1>
775{
777public:
785 template<class P = Policy>
786 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
787 {
788 return Base::intersection(geo2, geo1, intersection);
789 }
790};
791
796template <class Geometry1, class Geometry2, class Policy>
797class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 2>
798{
799 enum { dimworld = 3 };
800 enum { dim1 = 3 };
801 enum { dim2 = 2 };
802
803public:
804 using ctype = typename Policy::ctype;
805 using Point = typename Policy::Point;
806 using Intersection = typename Policy::Intersection;
807
808private:
809 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
810 using ReferenceElementsGeo1 = typename Dune::ReferenceElements<ctype, dim1>;
811 using ReferenceElementsGeo2 = typename Dune::ReferenceElements<ctype, dim2>;
812
813public:
828 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
829 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
830 {
831 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
832
833 // the candidate intersection points
834 std::vector<Point> points; points.reserve(10);
835
836 // add 3d geometry corners that are inside the 2d geometry
837 for (int i = 0; i < geo1.corners(); ++i)
838 if (intersectsPointGeometry(geo1.corner(i), geo2))
839 points.emplace_back(geo1.corner(i));
840
841 // add 2d geometry corners that are inside the 3d geometry
842 for (int i = 0; i < geo2.corners(); ++i)
843 if (intersectsPointGeometry(geo2.corner(i), geo1))
844 points.emplace_back(geo2.corner(i));
845
846 // get some geometry types
847 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
848 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
849
850 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
851 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
852
853 // intersection policy for point-like intersections (used later)
855
856 // add intersection points of all polyhedron edges (codim dim-1) with the polygon
857 for (int i = 0; i < referenceElement1.size(dim1-1); ++i)
858 {
859 const auto localEdgeGeom = referenceElement1.template geometry<dim1-1>(i);
860 const auto p = geo1.global(localEdgeGeom.corner(0));
861 const auto q = geo1.global(localEdgeGeom.corner(1));
862 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
863
865 typename PolySegTest::Intersection polySegIntersection;
866 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
867 points.emplace_back(polySegIntersection);
868 }
869
870 // add intersection points of all polygon faces (codim 1) with the polyhedron faces
871 for (int i = 0; i < referenceElement1.size(1); ++i)
872 {
873 const auto faceGeo = [&]()
874 {
875 const auto localFaceGeo = referenceElement1.template geometry<1>(i);
876 if (localFaceGeo.corners() == 4)
877 {
878 const auto a = geo1.global(localFaceGeo.corner(0));
879 const auto b = geo1.global(localFaceGeo.corner(1));
880 const auto c = geo1.global(localFaceGeo.corner(2));
881 const auto d = geo1.global(localFaceGeo.corner(3));
882
883 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
884 }
885 else
886 {
887 const auto a = geo1.global(localFaceGeo.corner(0));
888 const auto b = geo1.global(localFaceGeo.corner(1));
889 const auto c = geo1.global(localFaceGeo.corner(2));
890
891 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
892 }
893 }();
894
895 for (int j = 0; j < referenceElement2.size(1); ++j)
896 {
897 const auto localEdgeGeom = referenceElement2.template geometry<1>(j);
898 const auto p = geo2.global(localEdgeGeom.corner(0));
899 const auto q = geo2.global(localEdgeGeom.corner(1));
900
901 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
902
904 typename PolySegTest::Intersection polySegIntersection;
905 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
906 points.emplace_back(polySegIntersection);
907 }
908 }
909
910 // return if no intersection points were found
911 if (points.empty()) return false;
912
913 // remove duplicates
914 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
915 std::sort(points.begin(), points.end(), [&eps](const auto& a, const auto& b) -> bool
916 {
917 using std::abs;
918 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])));
919 });
920
921 auto removeIt = std::unique(points.begin(), points.end(), [&eps](const auto& a, const auto&b)
922 {
923 return (b-a).two_norm() < eps;
924 });
925
926 points.erase(removeIt, points.end());
927
928 // return false if we don't have more than three unique points
929 if (points.size() < 3) return false;
930
931 // intersection polygon is convex hull of above points
932 intersection = grahamConvexHull<2>(points);
933 assert(!intersection.empty());
934
935 return true;
936 }
937
952 template<class P = Policy, std::enable_if_t<P::dimIntersection != 2, int> = 0>
953 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
954 {
955 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
956 DUNE_THROW(Dune::NotImplemented, "Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
957 }
958};
959
964template <class Geometry1, class Geometry2, class Policy>
965class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 3>
966: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 2>
967{
969public:
977 template<class P = Policy>
978 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
979 {
980 return Base::intersection(geo2, geo1, intersection);
981 }
982};
983
988template <class Geometry1, class Geometry2, class Policy>
989class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 1>
990{
991 enum { dimworld = 3 };
992 enum { dim1 = 2 };
993 enum { dim2 = 1 };
994
995public:
996 using ctype = typename Policy::ctype;
997 using Point = typename Policy::Point;
998 using Intersection = typename Policy::Intersection;
999
1000private:
1001 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1002 using ReferenceElements = typename Dune::ReferenceElements<ctype, dim1>;
1003
1004public:
1014 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, 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
1019 ctype tfirst, tlast;
1020 if (intersect_(geo1, geo2, tfirst, tlast))
1021 {
1022 // the intersection exists. Export the intersections geometry now:
1023 // s(t) = a + t(b-a) in [tfirst, tlast]
1024 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
1025 return true;
1026 }
1027
1028 return false;
1029 }
1030
1040 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1041 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& is)
1042 {
1043 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1044
1045 const auto p = geo2.corner(0);
1046 const auto q = geo2.corner(1);
1047
1048 const auto a = geo1.corner(0);
1049 const auto b = geo1.corner(1);
1050 const auto c = geo1.corner(2);
1051
1052 if (geo1.corners() == 3)
1053 return intersect_<Policy>(a, b, c, p, q, is);
1054
1055 else if (geo1.corners() == 4)
1056 {
1057 Intersection is1, is2;
1058 bool hasSegment1, hasSegment2;
1059
1060 const auto d = geo1.corner(3);
1061 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1062 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1063
1064 if (hasSegment1 || hasSegment2)
1065 return false;
1066
1067 if (intersects1) { is = is1; return true; }
1068 if (intersects2) { is = is2; return true; }
1069
1070 return false;
1071 }
1072
1073 else
1074 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1075 << geo1.type() << ", "<< geo1.corners() << " corners.");
1076 }
1077
1078private:
1082 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1083 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
1084 {
1085 // lambda to obtain the facet corners on geo1
1086 auto getFacetCorners = [] (const Geometry1& geo1)
1087 {
1088 std::vector< std::array<int, 2> > facetCorners;
1089 switch (geo1.corners())
1090 {
1091 case 4: // quadrilateral
1092 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1093 break;
1094 case 3: // triangle
1095 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1096 break;
1097 default:
1098 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1099 << geo1.type() << " with "<< geo1.corners() << " corners.");
1100 }
1101
1102 return facetCorners;
1103 };
1104
1105 const auto center1 = geo1.center();
1106 const auto normal1 = crossProduct(geo1.corner(1) - geo1.corner(0),
1107 geo1.corner(2) - geo1.corner(0));
1108
1109 // lambda to obtain the normal on a facet
1110 auto computeNormal = [&center1, &normal1, &geo1] (const std::array<int, 2>& facetCorners)
1111 {
1112 const auto c0 = geo1.corner(facetCorners[0]);
1113 const auto c1 = geo1.corner(facetCorners[1]);
1114 const auto edge = c1 - c0;
1115
1116 Dune::FieldVector<ctype, dimworld> n = crossProduct(edge, normal1);
1117 n /= n.two_norm();
1118
1119 // orientation of the element is not known
1120 // make sure normal points outwards of element
1121 if ( n*(center1-c0) > 0.0 )
1122 n *= -1.0;
1123
1124 return n;
1125 };
1126
1127 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
1128 }
1129
1133 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1134 static bool intersect_(const Point& a, const Point& b, const Point& c,
1135 const Point& p, const Point& q,
1136 Intersection& is)
1137 {
1138 bool hasSegment;
1139 return intersect_(a, b, c, p, q, is, hasSegment);
1140 }
1141
1146 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1147 static bool intersect_(const Point& a, const Point& b, const Point& c,
1148 const Point& p, const Point& q,
1149 Intersection& is, bool& hasSegmentIs)
1150 {
1151 hasSegmentIs = false;
1152
1153 const auto ab = b - a;
1154 const auto ac = c - a;
1155 const auto qp = p - q;
1156
1157 // compute the triangle normal that defines the triangle plane
1158 const auto normal = crossProduct(ab, ac);
1159
1160 // compute the denominator
1161 // if denom is 0 the segment is parallel and we can return
1162 const auto denom = normal*qp;
1163 const auto abnorm2 = ab.two_norm2();
1164 const auto eps = eps_*abnorm2*qp.two_norm();
1165 using std::abs;
1166 if (abs(denom) < eps)
1167 {
1168 const auto pa = a - p;
1169 const auto denom2 = normal*pa;
1170 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1171 return false;
1172
1173 // if we get here, we are in-plane. Check if a
1174 // segment-like intersection with geometry 1 exists.
1175 using SegmentPolicy = typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1176 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1177 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1178 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1179 using SegmentIntersectionType = typename SegmentIntersectionAlgorithm::Intersection;
1180 SegmentIntersectionType segmentIs;
1181
1182 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1183 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1184 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1185 {
1186 hasSegmentIs = true;
1187 return false;
1188 }
1189
1190 // We are in-plane. A point-like
1191 // intersection can only be on the edges
1192 for (const auto& ip : {p, q})
1193 {
1194 if ( intersectsPointSimplex(ip, a, b)
1195 || intersectsPointSimplex(ip, b, c)
1196 || intersectsPointSimplex(ip, c, a) )
1197 {
1198 is = ip;
1199 return true;
1200 }
1201 }
1202
1203 return false;
1204 }
1205
1206 // compute intersection t value of pq with plane of triangle.
1207 // a segment intersects if and only if 0 <= t <= 1.
1208 const auto ap = p - a;
1209 const auto t = (ap*normal)/denom;
1210 if (t < 0.0 - eps_) return false;
1211 if (t > 1.0 + eps_) return false;
1212
1213 // compute the barycentric coordinates and check if the intersection point
1214 // is inside the bounds of the triangle
1215 const auto e = crossProduct(qp, ap);
1216 const auto v = (ac*e)/denom;
1217 if (v < -eps_ || v > 1.0 + eps_) return false;
1218 const auto w = -(ab*e)/denom;
1219 if (w < -eps_ || v + w > 1.0 + eps_) return false;
1220
1221 // Now we are sure there is an intersection points
1222 // Perform delayed division compute the last barycentric coordinate component
1223 const auto u = 1.0 - v - w;
1224
1225 Point ip(0.0);
1226 ip.axpy(u, a);
1227 ip.axpy(v, b);
1228 ip.axpy(w, c);
1229 is = ip;
1230 return true;
1231 }
1232};
1233
1238template <class Geometry1, class Geometry2, class Policy>
1239class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 2>
1240: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 2, 1>
1241{
1243public:
1251 template<class P = Policy>
1252 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
1253 {
1254 return Base::intersection(geo2, geo1, intersection);
1255 }
1256};
1257
1262template <class Geometry1, class Geometry2, class Policy>
1263class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 1>
1264{
1265 enum { dimworld = 3 };
1266 enum { dim1 = 1 };
1267 enum { dim2 = 1 };
1268
1269public:
1270 using ctype = typename Policy::ctype;
1271 using Point = typename Policy::Point;
1272 using Intersection = typename Policy::Intersection;
1273
1274private:
1275 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1276 using ReferenceElements = typename Dune::ReferenceElements<ctype, dim1>;
1277
1278public:
1286 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1287 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1288 {
1289 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1290 DUNE_THROW(Dune::NotImplemented, "segment-segment intersection detection for point-like intersections");
1291 }
1292
1300 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1301 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1302 {
1303 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1304
1305 const auto& a = geo1.corner(0);
1306 const auto& b = geo1.corner(1);
1307 const auto& p = geo2.corner(0);
1308 const auto& q = geo2.corner(1);
1309
1310 const auto ab = b-a;
1311 const auto pq = q-p;
1312 const auto abNorm = ab.two_norm();
1313 const auto pqNorm = pq.two_norm();
1314
1315 using std::max;
1316 const auto maxNorm = max(abNorm, pqNorm);
1317 const auto eps2 = eps_*maxNorm*maxNorm;
1318
1319 // if the segment are not parallel there is no segment intersection
1320 using std::abs;
1321 if (!(abs(abs(ab*pq) - abNorm*pqNorm) < eps2))
1322 return false;
1323
1324 const auto ap = (p-a);
1325 const auto aq = (q-a);
1326
1327 // points have to be colinear
1328 if (!(crossProduct(ap, aq).two_norm() < eps2))
1329 return false;
1330
1331 // scaled local coordinates
1332 // we save the division for the normalization for now
1333 // and do it in the very end if we find an intersection
1334 auto tp = ap*ab;
1335 auto tq = aq*ab;
1336 const auto abnorm2 = abNorm*abNorm;
1337
1338 // make sure they are sorted
1339 using std::swap;
1340 if (tp > tq)
1341 swap(tp, tq);
1342
1343 using std::min; using std::max;
1344 tp = min(abnorm2, max(0.0, tp));
1345 tq = max(0.0, min(abnorm2, tq));
1346
1347 if (abs(tp-tq) < eps2)
1348 return false;
1349
1350 intersection = Intersection({geo1.global(tp/abnorm2), geo1.global(tq/abnorm2)});
1351 return true;
1352 }
1353};
1354
1355} // end namespace Dumux
1356
1357# endif
An axis-aligned bounding box volume hierarchy for dune grids.
Detect if a point intersects a geometry.
A function to compute the convex hull of a point cloud and a function to triangulate the polygon area...
Define some often used mathematical functions.
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:631
bool intersectsPointGeometry(const Dune::FieldVector< ctype, dimworld > &point, const Geometry &g)
Find out whether a point is inside a three-dimensional geometry.
Definition: intersectspointgeometry.hh:38
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: intersectspointsimplex.hh:36
Definition: adapt.hh:29
bool intersectsPointBoundingBox(const Dune::FieldVector< ctype, dimworld > &point, const ctype *b)
Check whether a point is intersectin a bounding box (dimworld == 3)
Definition: boundingboxtree.hh:304
typename DefaultPolicyChooser< Geometry1, Geometry2 >::type DefaultPolicy
Helper alias to define the default intersection policy.
Definition: geometryintersection.hh:114
Policy structure for point-like intersections.
Definition: geometryintersection.hh:44
static constexpr int dimIntersection
Definition: geometryintersection.hh:49
Point Intersection
Definition: geometryintersection.hh:50
ct ctype
Definition: geometryintersection.hh:45
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:46
static constexpr int dimWorld
Definition: geometryintersection.hh:48
Policy structure for segment-like intersections.
Definition: geometryintersection.hh:56
Segment Intersection
Definition: geometryintersection.hh:63
static constexpr int dimIntersection
Definition: geometryintersection.hh:62
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:58
static constexpr int dimWorld
Definition: geometryintersection.hh:61
std::array< Point, 2 > Segment
Definition: geometryintersection.hh:59
ct ctype
Definition: geometryintersection.hh:57
Policy structure for polygon-like intersections.
Definition: geometryintersection.hh:69
static constexpr int dimWorld
Definition: geometryintersection.hh:74
static constexpr int dimIntersection
Definition: geometryintersection.hh:75
ct ctype
Definition: geometryintersection.hh:70
Polygon Intersection
Definition: geometryintersection.hh:76
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:71
std::vector< Point > Polygon
Definition: geometryintersection.hh:72
Policy structure for polyhedron-like intersections.
Definition: geometryintersection.hh:82
std::vector< Point > Polyhedron
Definition: geometryintersection.hh:85
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:84
Polyhedron Intersection
Definition: geometryintersection.hh:89
static constexpr int dimIntersection
Definition: geometryintersection.hh:88
static constexpr int dimWorld
Definition: geometryintersection.hh:87
ct ctype
Definition: geometryintersection.hh:83
default policy chooser class
Definition: geometryintersection.hh:95
std::tuple_element_t< isDim, DefaultPolicies > type
Definition: geometryintersection.hh:109
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:217
typename Policy::ctype ctype
Definition: geometryintersection.hh:221
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:223
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:226
typename Policy::Point Point
Definition: geometryintersection.hh:222
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:251
typename Policy::ctype ctype
Definition: geometryintersection.hh:249
typename Policy::Point Point
Definition: geometryintersection.hh:250
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:260
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:380
typename Policy::ctype ctype
Definition: geometryintersection.hh:363
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:365
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:404
typename Policy::Point Point
Definition: geometryintersection.hh:364
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:491
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:511
typename Policy::ctype ctype
Definition: geometryintersection.hh:509
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:532
typename Policy::Point Point
Definition: geometryintersection.hh:510
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:657
typename Policy::Point Point
Definition: geometryintersection.hh:656
typename Policy::ctype ctype
Definition: geometryintersection.hh:655
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:676
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:786
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:806
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:829
typename Policy::Point Point
Definition: geometryintersection.hh:805
typename Policy::ctype ctype
Definition: geometryintersection.hh:804
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:978
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:998
typename Policy::Point Point
Definition: geometryintersection.hh:997
typename Policy::ctype ctype
Definition: geometryintersection.hh:996
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1015
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1041
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:1252
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:1272
typename Policy::Point Point
Definition: geometryintersection.hh:1271
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:1287
typename Policy::ctype ctype
Definition: geometryintersection.hh:1270