3.1-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 <dune/common/exceptions.hh>
27#include <dune/common/promotiontraits.hh>
28#include <dune/geometry/referenceelements.hh>
29#include <dune/geometry/multilineargeometry.hh>
30
31#include <dumux/common/math.hh>
35
36namespace Dumux {
37namespace IntersectionPolicy {
38
40template<class ct, int dw>
42{
43 using ctype = ct;
44 using Point = Dune::FieldVector<ctype, dw>;
45
46 static constexpr int dimWorld = dw;
47 static constexpr int dimIntersection = 0;
49};
50
52template<class ct, int dw>
54{
55 using ctype = ct;
56 using Point = Dune::FieldVector<ctype, dw>;
57 using Segment = std::array<Point, 2>;
58
59 static constexpr int dimWorld = dw;
60 static constexpr int dimIntersection = 1;
62};
63
65template<class ct, int dw>
67{
68 using ctype = ct;
69 using Point = Dune::FieldVector<ctype, dw>;
70 using Polygon = std::vector<Point>;
71
72 static constexpr int dimWorld = dw;
73 static constexpr int dimIntersection = 2;
75};
76
78template<class Geometry1, class Geometry2, int dim2 = Geometry2::mydimension>
80{
81 static constexpr int dimworld = Geometry1::coorddimension;
82 static constexpr int isDim = std::min( int(Geometry1::mydimension), int(Geometry2::mydimension) );
83 static_assert(dimworld == int(Geometry2::coorddimension), "Geometries must have the same coordinate dimension!");
84 static_assert(isDim == 1 || isDim == 2, "No chooser class specialization available for given dimensions");
85 using ctype = typename Dune::PromotionTraits<typename Geometry1::ctype, typename Geometry2::ctype>::PromotedType;
86public:
87 using type = std::conditional_t< isDim == 2,
90};
91
93template<class Geometry1, class Geometry2>
95
96} // end namespace IntersectionPolicy
97
98namespace Impl {
99
113template< class Geo1, class Geo2, class ctype,
114 class GetFacetCornerIndices, class ComputeNormalFunction >
115bool computeSegmentIntersection(const Geo1& geo1, const Geo2& geo2, ctype baseEps,
116 ctype& tfirst, ctype& tlast,
117 const GetFacetCornerIndices& getFacetCornerIndices,
118 const ComputeNormalFunction& computeNormal)
119{
120 static_assert(int(Geo2::mydimension) == 1, "Geometry2 must be a segment");
121 static_assert(int(Geo1::mydimension) > int(Geo2::mydimension),
122 "Dimension of Geometry1 must be higher than that of Geometry2");
123
124 const auto a = geo2.corner(0);
125 const auto b = geo2.corner(1);
126 const auto d = b - a;
127
128 // The initial interval is the whole segment.
129 // Afterwards we start clipping the interval
130 // at the intersections with the facets of geo1
131 tfirst = 0.0;
132 tlast = 1.0;
133
134 const auto facets = getFacetCornerIndices(geo1);
135 for (const auto& f : facets)
136 {
137 const auto n = computeNormal(f);
138
139 const auto c0 = geo1.corner(f[0]);
140 const ctype denom = n*d;
141 const ctype dist = n*(a-c0);
142
143 // use first edge of the facet to scale eps
144 const auto edge1 = geo1.corner(f[1]) - geo1.corner(f[0]);
145 const ctype eps = baseEps*edge1.two_norm();
146
147 // if denominator is zero the segment in parallel to the edge.
148 // If the distance is positive there is no intersection
149 using std::abs;
150 if (abs(denom) < eps)
151 {
152 if (dist > eps)
153 return false;
154 }
155 else // not parallel: compute line-line intersection
156 {
157 const ctype t = -dist / denom;
158 // if entering half space cut tfirst if t is larger
159 using std::signbit;
160 if (signbit(denom))
161 {
162 if (t > tfirst)
163 tfirst = t;
164 }
165 // if exiting half space cut tlast if t is smaller
166 else
167 {
168 if (t < tlast)
169 tlast = t;
170 }
171 // there is no intersection of the interval is empty
172 // use unscaled epsilon since t is in local coordinates
173 if (tfirst > tlast - baseEps)
174 return false;
175 }
176 }
177
178 // If we made it until here an intersection exists.
179 return true;
180}
181
182} // end namespace Impl
183
190template
191<class Geometry1, class Geometry2,
192 class Policy = IntersectionPolicy::DefaultPolicy<Geometry1, Geometry2>,
193 int dimworld = Geometry1::coorddimension,
194 int dim1 = Geometry1::mydimension,
195 int dim2 = Geometry2::mydimension>
197{
198 static constexpr int dimWorld = Policy::dimWorld;
199
200public:
201 using ctype = typename Policy::ctype;
202 using Point = typename Policy::Point;
203 using Intersection = typename Policy::Intersection;
204
206 using IntersectionType [[deprecated("Please use Intersection instead")]] = Intersection;
207
209 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
210 {
211 static_assert(dimworld == Geometry2::coorddimension, "Can only intersect geometries of same coordinate dimension");
212 DUNE_THROW(Dune::NotImplemented, "Geometry intersection detection with intersection computation not implemented for dimworld = "
213 << dimworld << ", dim1 = " << dim1 << ", dim2 = " << dim2);
214 }
215};
216
221template <class Geometry1, class Geometry2, class Policy>
222class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 1>
223{
224 enum { dimworld = 2 };
225 enum { dim1 = 1 };
226 enum { dim2 = 1 };
227
228 // base epsilon for floating point comparisons
229 static constexpr typename Policy::ctype eps_ = 1.5e-7;
230
231public:
232 using ctype = typename Policy::ctype;
233 using Point = typename Policy::Point;
234 using Intersection = typename Policy::Intersection;
235
237 using IntersectionType [[deprecated("Please use Intersection instead")]] = Intersection;
238
245 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
246 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
247 {
248 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
249
250 const auto v1 = geo1.corner(1) - geo1.corner(0);
251 const auto v2 = geo2.corner(1) - geo2.corner(0);
252 const auto ac = geo2.corner(0) - geo1.corner(0);
253
254 auto n2 = Point({-1.0*v2[1], v2[0]});
255 n2 /= n2.two_norm();
256
257 // compute distance of first corner on geo2 to line1
258 const auto dist12 = n2*ac;
259
260 // first check parallel segments
261 using std::abs;
262 using std::sqrt;
263
264 const auto v1Norm2 = v1.two_norm2();
265 const auto eps = eps_*sqrt(v1Norm2);
266 const auto eps2 = eps_*v1Norm2;
267
268 const auto sp = n2*v1;
269 if (abs(sp) < eps)
270 {
271 if (abs(dist12) > eps)
272 return false;
273
274 // point intersection can only be one of the corners
275 if ( (v1*v2) < 0.0 )
276 {
277 if ( ac.two_norm2() < eps2 )
278 { intersection = geo2.corner(0); return true; }
279
280 if ( (geo2.corner(1) - geo1.corner(1)).two_norm2() < eps2 )
281 { intersection = geo2.corner(1); return true; }
282 }
283 else
284 {
285 if ( (geo2.corner(1) - geo1.corner(0)).two_norm2() < eps2 )
286 { intersection = geo2.corner(1); return true; }
287
288 if ( (geo2.corner(0) - geo1.corner(1)).two_norm2() < eps2 )
289 { intersection = geo2.corner(0); return true; }
290 }
291
292 // no intersection
293 return false;
294 }
295
296 // intersection point on v1 in local coords
297 const auto t1 = dist12 / sp;
298
299 // check if the local coords are valid
300 if (t1 < -1.0*eps_ || t1 > 1.0 + eps_)
301 return false;
302
303 // compute global coordinates
304 auto isPoint = geo1.global(t1);
305
306 // check if point is in bounding box of v2
307 const auto c = geo2.corner(0);
308 const auto d = geo2.corner(1);
309
310 using std::min; using std::max;
311 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]) });
312
313 if ( intersectsPointBoundingBox(isPoint, bBox.data()) )
314 {
315 intersection = std::move(isPoint);
316 return true;
317 }
318
319 return false;
320 }
321
329 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
330 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
331 {
332 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
333 DUNE_THROW(Dune::NotImplemented, "segment-segment intersection detection for segment-like intersections");
334 }
335};
336
341template <class Geometry1, class Geometry2, class Policy>
342class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 1>
343{
344 enum { dimworld = 2 };
345 enum { dim1 = 2 };
346 enum { dim2 = 1 };
347
348public:
349 using ctype = typename Policy::ctype;
350 using Point = typename Policy::Point;
351 using Intersection = typename Policy::Intersection;
352
354 using IntersectionType [[deprecated("Please use Intersection instead")]] = Intersection;
355
356private:
357 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
358 using ReferenceElements = typename Dune::ReferenceElements<ctype, dim1>;
359
360public:
368 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
369 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
370 {
371 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
372
373 ctype tfirst, tlast;
374 if (intersect_(geo1, geo2, tfirst, tlast))
375 {
376 // the intersection exists. Export the intersections geometry now:
377 // s(t) = a + t(b-a) in [tfirst, tlast]
378 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
379 return true;
380 }
381
382 return false;
383 }
384
392 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
393 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename P::Intersection& intersection)
394 {
395 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
396
397 ctype tfirst, tlast;
398 if (!intersect_(geo1, geo2, tfirst, tlast))
399 {
400 // if start & end point are same, the intersection is a point
401 using std::abs;
402 if ( tfirst > tlast - eps_ )
403 {
404 intersection = Intersection(geo2.global(tfirst));
405 return true;
406 }
407 }
408
409 return false;
410 }
411
412private:
416 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
417 {
418 // lambda to obtain the facet corners on geo1
419 auto getFacetCorners = [] (const Geometry1& geo1)
420 {
421 std::vector< std::array<int, 2> > facetCorners;
422 switch (geo1.corners())
423 {
424 case 4: // quadrilateral
425 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
426 break;
427 case 3: // triangle
428 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
429 break;
430 default:
431 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
432 << geo1.type() << " with "<< geo1.corners() << " corners.");
433 }
434
435 return facetCorners;
436 };
437
438 // lambda to obtain the normal on a facet
439 const auto center1 = geo1.center();
440 auto computeNormal = [&geo1, &center1] (const std::array<int, 2>& facetCorners)
441 {
442 const auto c0 = geo1.corner(facetCorners[0]);
443 const auto c1 = geo1.corner(facetCorners[1]);
444 const auto edge = c1 - c0;
445
446 Dune::FieldVector<ctype, dimworld> n({-1.0*edge[1], edge[0]});
447 n /= n.two_norm();
448
449 // orientation of the element is not known
450 // make sure normal points outwards of element
451 if ( n*(center1-c0) > 0.0 )
452 n *= -1.0;
453
454 return n;
455 };
456
457 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
458 }
459};
460
465template <class Geometry1, class Geometry2, class Policy>
466class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 1, 2>
467: public GeometryIntersection<Geometry2, Geometry1, Policy, 2, 2, 1>
468{
470
471public:
479 template<class P = Policy>
480 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
481 {
482 return Base::intersection(geo2, geo1, intersection);
483 }
484};
485
490template <class Geometry1, class Geometry2, class Policy>
491class GeometryIntersection<Geometry1, Geometry2, Policy, 2, 2, 2>
492{
493 enum { dimworld = 2 };
494 enum { dim1 = 2 };
495 enum { dim2 = 2 };
496
497public:
498 using ctype = typename Policy::ctype;
499 using Point = typename Policy::Point;
500 using Intersection = typename Policy::Intersection;
501
503 using IntersectionType [[deprecated("Please use Intersection instead")]] = Intersection;
504
505private:
506 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
507 using ReferenceElementsGeo1 = typename Dune::ReferenceElements<ctype, dim1>;
508 using ReferenceElementsGeo2 = typename Dune::ReferenceElements<ctype, dim2>;
509
510public:
523 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
524 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
525 {
526 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
527
528 // the candidate intersection points
529 std::vector<Point> points; points.reserve(6);
530
531 // add polygon1 corners that are inside polygon2
532 for (int i = 0; i < geo1.corners(); ++i)
533 if (intersectsPointGeometry(geo1.corner(i), geo2))
534 points.emplace_back(geo1.corner(i));
535
536 const auto numPoints1 = points.size();
537 if (numPoints1 != geo1.corners())
538 {
539 // add polygon2 corners that are inside polygon1
540 for (int i = 0; i < geo2.corners(); ++i)
541 if (intersectsPointGeometry(geo2.corner(i), geo1))
542 points.emplace_back(geo2.corner(i));
543
544 if (points.empty())
545 return false;
546
547 if (points.size() - numPoints1 != geo2.corners())
548 {
549 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
550 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
551
552 // add intersections of edges
553 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
555 for (int i = 0; i < referenceElement1.size(dim1-1); ++i)
556 {
557 const auto localEdgeGeom1 = referenceElement1.template geometry<dim1-1>(i);
558 const auto edge1 = SegGeometry( Dune::GeometryTypes::line,
559 std::vector<Point>( {geo1.global(localEdgeGeom1.corner(0)),
560 geo1.global(localEdgeGeom1.corner(1))} ));
561
562 for (int j = 0; j < referenceElement2.size(dim2-1); ++j)
563 {
564 const auto localEdgeGeom2 = referenceElement2.template geometry<dim2-1>(j);
565 const auto edge2 = SegGeometry( Dune::GeometryTypes::line,
566 std::vector<Point>( {geo2.global(localEdgeGeom2.corner(0)),
567 geo2.global(localEdgeGeom2.corner(1))} ));
568
570 typename EdgeTest::Intersection edgeIntersection;
571 if (EdgeTest::intersection(edge1, edge2, edgeIntersection))
572 points.emplace_back(edgeIntersection);
573 }
574 }
575 }
576 }
577
578 if (points.empty())
579 return false;
580
581 // remove duplicates
582 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
583 std::sort(points.begin(), points.end(), [&eps](const auto& a, const auto& b) -> bool
584 {
585 using std::abs;
586 return (abs(a[0]-b[0]) > eps ? a[0] < b[0] : a[1] < b[1]);
587 });
588
589 auto removeIt = std::unique(points.begin(), points.end(), [&eps](const auto& a, const auto&b)
590 {
591 return (b-a).two_norm() < eps;
592 });
593
594 points.erase(removeIt, points.end());
595
596 // return false if we don't have at least three unique points
597 if (points.size() < 3)
598 return false;
599
600 // intersection polygon is convex hull of above points
601 intersection = grahamConvexHull<2>(points);
602 assert(!intersection.empty());
603 return true;
604 }
605
613 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
614 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
615 {
616 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
617 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for segment-like intersections");
618 }
619
627 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
628 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
629 {
630 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
631 DUNE_THROW(Dune::NotImplemented, "Polygon-polygon intersection detection for touching points");
632 }
633};
634
639template <class Geometry1, class Geometry2, class Policy>
640class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 1>
641{
642 enum { dimworld = 3 };
643 enum { dim1 = 3 };
644 enum { dim2 = 1 };
645
646public:
647 using ctype = typename Policy::ctype;
648 using Point = typename Policy::Point;
649 using Intersection = typename Policy::Intersection;
650
652 using IntersectionType [[deprecated("Please use Intersection instead")]] = Intersection;
653
654private:
655 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
656 using ReferenceElements = typename Dune::ReferenceElements<ctype, dim1>;
657
658public:
670 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
671 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
672 {
673 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
674
675 ctype tfirst, tlast;
676 if (intersect_(geo1, geo2, tfirst, tlast))
677 {
678 // Intersection exists. Export the intersections geometry now:
679 // s(t) = a + t(b-a) in [tfirst, tlast]
680 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
681 return true;
682 }
683
684 return false;
685 }
686
698 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
699 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
700 {
701 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
702
703 ctype tfirst, tlast;
704 if (intersect_(geo1, geo2, tfirst, tlast))
705 {
706 // if start & end point are same, the intersection is a point
707 using std::abs;
708 if ( abs(tfirst - tlast) < eps_ )
709 {
710 intersection = Intersection(geo2.global(tfirst));
711 return true;
712 }
713 }
714
715 return false;
716 }
717
718private:
722 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
723 {
724 // lambda to obtain facet corners on geo1
725 auto getFacetCorners = [] (const Geometry1& geo1)
726 {
727 std::vector< std::vector<int> > facetCorners;
728 // sort facet corner so that normal n = (p1-p0)x(p2-p0) always points outwards
729 switch (geo1.corners())
730 {
731 // todo compute intersection geometries
732 case 8: // hexahedron
733 facetCorners = {{2, 0, 6, 4}, {1, 3, 5, 7}, {0, 1, 4, 5},
734 {3, 2, 7, 6}, {1, 0, 3, 2}, {4, 5, 6, 7}};
735 break;
736 case 4: // tetrahedron
737 facetCorners = {{1, 0, 2}, {0, 1, 3}, {0, 3, 2}, {1, 2, 3}};
738 break;
739 default:
740 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
741 << geo1.type() << ", "<< geo1.corners() << " corners.");
742 }
743
744 return facetCorners;
745 };
746
747 // lambda to obtain the normal on a facet
748 auto computeNormal = [&geo1] (const std::vector<int>& facetCorners)
749 {
750 const auto v0 = geo1.corner(facetCorners[1]) - geo1.corner(facetCorners[0]);
751 const auto v1 = geo1.corner(facetCorners[2]) - geo1.corner(facetCorners[0]);
752
753 auto n = crossProduct(v0, v1);
754 n /= n.two_norm();
755
756 return n;
757 };
758
759 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
760 }
761};
762
767template <class Geometry1, class Geometry2, class Policy>
768class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 3>
769: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 1>
770{
772public:
780 template<class P = Policy>
781 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
782 {
783 return Base::intersection(geo2, geo1, intersection);
784 }
785};
786
791template <class Geometry1, class Geometry2, class Policy>
792class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 3, 2>
793{
794 enum { dimworld = 3 };
795 enum { dim1 = 3 };
796 enum { dim2 = 2 };
797
798public:
799 using ctype = typename Policy::ctype;
800 using Point = typename Policy::Point;
801 using Intersection = typename Policy::Intersection;
802
804 using IntersectionType [[deprecated("Please use Intersection instead")]] = Intersection;
805
806private:
807 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
808 using ReferenceElementsGeo1 = typename Dune::ReferenceElements<ctype, dim1>;
809 using ReferenceElementsGeo2 = typename Dune::ReferenceElements<ctype, dim2>;
810
811public:
826 template<class P = Policy, std::enable_if_t<P::dimIntersection == 2, int> = 0>
827 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
828 {
829 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
830
831 // the candidate intersection points
832 std::vector<Point> points; points.reserve(10);
833
834 // add 3d geometry corners that are inside the 2d geometry
835 for (int i = 0; i < geo1.corners(); ++i)
836 if (intersectsPointGeometry(geo1.corner(i), geo2))
837 points.emplace_back(geo1.corner(i));
838
839 // add 2d geometry corners that are inside the 3d geometry
840 for (int i = 0; i < geo2.corners(); ++i)
841 if (intersectsPointGeometry(geo2.corner(i), geo1))
842 points.emplace_back(geo2.corner(i));
843
844 // get some geometry types
845 using PolyhedronFaceGeometry = Dune::MultiLinearGeometry<ctype, 2, dimworld>;
846 using SegGeometry = Dune::MultiLinearGeometry<ctype, 1, dimworld>;
847
848 const auto referenceElement1 = ReferenceElementsGeo1::general(geo1.type());
849 const auto referenceElement2 = ReferenceElementsGeo2::general(geo2.type());
850
851 // intersection policy for point-like intersections (used later)
853
854 // add intersection points of all polyhedron edges (codim dim-1) with the polygon
855 for (int i = 0; i < referenceElement1.size(dim1-1); ++i)
856 {
857 const auto localEdgeGeom = referenceElement1.template geometry<dim1-1>(i);
858 const auto p = geo1.global(localEdgeGeom.corner(0));
859 const auto q = geo1.global(localEdgeGeom.corner(1));
860 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
861
863 typename PolySegTest::Intersection polySegIntersection;
864 if (PolySegTest::intersection(geo2, segGeo, polySegIntersection))
865 points.emplace_back(polySegIntersection);
866 }
867
868 // add intersection points of all polygon faces (codim 1) with the polyhedron faces
869 for (int i = 0; i < referenceElement1.size(1); ++i)
870 {
871 const auto faceGeo = [&]()
872 {
873 const auto localFaceGeo = referenceElement1.template geometry<1>(i);
874 if (localFaceGeo.corners() == 4)
875 {
876 const auto a = geo1.global(localFaceGeo.corner(0));
877 const auto b = geo1.global(localFaceGeo.corner(1));
878 const auto c = geo1.global(localFaceGeo.corner(2));
879 const auto d = geo1.global(localFaceGeo.corner(3));
880
881 return PolyhedronFaceGeometry(Dune::GeometryTypes::cube(2), std::vector<Point>{a, b, c, d});
882 }
883 else
884 {
885 const auto a = geo1.global(localFaceGeo.corner(0));
886 const auto b = geo1.global(localFaceGeo.corner(1));
887 const auto c = geo1.global(localFaceGeo.corner(2));
888
889 return PolyhedronFaceGeometry(Dune::GeometryTypes::simplex(2), std::vector<Point>{a, b, c});
890 }
891 }();
892
893 for (int j = 0; j < referenceElement2.size(1); ++j)
894 {
895 const auto localEdgeGeom = referenceElement2.template geometry<1>(j);
896 const auto p = geo2.global(localEdgeGeom.corner(0));
897 const auto q = geo2.global(localEdgeGeom.corner(1));
898
899 const auto segGeo = SegGeometry(Dune::GeometryTypes::line, std::vector<Point>{p, q});
900
902 typename PolySegTest::Intersection polySegIntersection;
903 if (PolySegTest::intersection(faceGeo, segGeo, polySegIntersection))
904 points.emplace_back(polySegIntersection);
905 }
906 }
907
908 // return if no intersection points were found
909 if (points.empty()) return false;
910
911 // remove duplicates
912 const auto eps = (geo1.corner(0) - geo1.corner(1)).two_norm()*eps_;
913 std::sort(points.begin(), points.end(), [&eps](const auto& a, const auto& b) -> bool
914 {
915 using std::abs;
916 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])));
917 });
918
919 auto removeIt = std::unique(points.begin(), points.end(), [&eps](const auto& a, const auto&b)
920 {
921 return (b-a).two_norm() < eps;
922 });
923
924 points.erase(removeIt, points.end());
925
926 // return false if we don't have more than three unique points
927 if (points.size() < 3) return false;
928
929 // intersection polygon is convex hull of above points
930 intersection = grahamConvexHull<2>(points);
931 assert(!intersection.empty());
932
933 return true;
934 }
935
950 template<class P = Policy, std::enable_if_t<P::dimIntersection != 2, int> = 0>
951 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
952 {
953 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
954 DUNE_THROW(Dune::NotImplemented, "Polyhedron-polygon intersection detection only implemented for polygon-like intersections");
955 }
956};
957
962template <class Geometry1, class Geometry2, class Policy>
963class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 3>
964: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 3, 2>
965{
967public:
975 template<class P = Policy>
976 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
977 {
978 return Base::intersection(geo2, geo1, intersection);
979 }
980};
981
986template <class Geometry1, class Geometry2, class Policy>
987class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 2, 1>
988{
989 enum { dimworld = 3 };
990 enum { dim1 = 2 };
991 enum { dim2 = 1 };
992
993public:
994 using ctype = typename Policy::ctype;
995 using Point = typename Policy::Point;
996 using Intersection = typename Policy::Intersection;
997
999 using IntersectionType [[deprecated("Please use Intersection instead")]] = Intersection;
1000
1001private:
1002 static constexpr ctype eps_ = 1.5e-7; // base epsilon for floating point comparisons
1003 using ReferenceElements = typename Dune::ReferenceElements<ctype, dim1>;
1004
1005public:
1015 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1016 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& intersection)
1017 {
1018 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1019
1020 ctype tfirst, tlast;
1021 if (intersect_(geo1, geo2, tfirst, tlast))
1022 {
1023 // the intersection exists. Export the intersections geometry now:
1024 // s(t) = a + t(b-a) in [tfirst, tlast]
1025 intersection = Intersection({geo2.global(tfirst), geo2.global(tlast)});
1026 return true;
1027 }
1028
1029 return false;
1030 }
1031
1041 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1042 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, Intersection& is)
1043 {
1044 static_assert(int(dimworld) == int(Geometry2::coorddimension), "Can only collide geometries of same coordinate dimension");
1045
1046 const auto p = geo2.corner(0);
1047 const auto q = geo2.corner(1);
1048
1049 const auto a = geo1.corner(0);
1050 const auto b = geo1.corner(1);
1051 const auto c = geo1.corner(2);
1052
1053 if (geo1.corners() == 3)
1054 return intersect_<Policy>(a, b, c, p, q, is);
1055
1056 else if (geo1.corners() == 4)
1057 {
1058 Intersection is1, is2;
1059 bool hasSegment1, hasSegment2;
1060
1061 const auto d = geo1.corner(3);
1062 const bool intersects1 = intersect_<Policy>(a, b, d, p, q, is1, hasSegment1);
1063 const bool intersects2 = intersect_<Policy>(a, d, c, p, q, is2, hasSegment2);
1064
1065 if (hasSegment1 || hasSegment2)
1066 return false;
1067
1068 if (intersects1) { is = is1; return true; }
1069 if (intersects2) { is = is2; return true; }
1070
1071 return false;
1072 }
1073
1074 else
1075 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1076 << geo1.type() << ", "<< geo1.corners() << " corners.");
1077 }
1078
1089 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1090 [[deprecated("Please use intersection(triangle, segment, ...) instead")]]
1091 static bool intersection(const Point& a, const Point& b, const Point& c,
1092 const Point& p, const Point& q,
1093 Intersection& is)
1094 {
1095 return intersect_<Policy>(a, b, c, p, q, is);
1096 }
1097
1098private:
1102 template<class P = Policy, std::enable_if_t<P::dimIntersection == 1, int> = 0>
1103 static bool intersect_(const Geometry1& geo1, const Geometry2& geo2, ctype& tfirst, ctype& tlast)
1104 {
1105 // lambda to obtain the facet corners on geo1
1106 auto getFacetCorners = [] (const Geometry1& geo1)
1107 {
1108 std::vector< std::array<int, 2> > facetCorners;
1109 switch (geo1.corners())
1110 {
1111 case 4: // quadrilateral
1112 facetCorners = {{0, 2}, {3, 1}, {1, 0}, {2, 3}};
1113 break;
1114 case 3: // triangle
1115 facetCorners = {{1, 0}, {0, 2}, {2, 1}};
1116 break;
1117 default:
1118 DUNE_THROW(Dune::NotImplemented, "Collision of segment and geometry of type "
1119 << geo1.type() << " with "<< geo1.corners() << " corners.");
1120 }
1121
1122 return facetCorners;
1123 };
1124
1125 const auto center1 = geo1.center();
1126 const auto normal1 = crossProduct(geo1.corner(1) - geo1.corner(0),
1127 geo1.corner(2) - geo1.corner(0));
1128
1129 // lambda to obtain the normal on a facet
1130 auto computeNormal = [&center1, &normal1, &geo1] (const std::array<int, 2>& facetCorners)
1131 {
1132 const auto c0 = geo1.corner(facetCorners[0]);
1133 const auto c1 = geo1.corner(facetCorners[1]);
1134 const auto edge = c1 - c0;
1135
1136 Dune::FieldVector<ctype, dimworld> n = crossProduct(edge, normal1);
1137 n /= n.two_norm();
1138
1139 // orientation of the element is not known
1140 // make sure normal points outwards of element
1141 if ( n*(center1-c0) > 0.0 )
1142 n *= -1.0;
1143
1144 return n;
1145 };
1146
1147 return Impl::computeSegmentIntersection(geo1, geo2, eps_, tfirst, tlast, getFacetCorners, computeNormal);
1148 }
1149
1153 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1154 static bool intersect_(const Point& a, const Point& b, const Point& c,
1155 const Point& p, const Point& q,
1156 Intersection& is)
1157 {
1158 bool hasSegment;
1159 return intersect_(a, b, c, p, q, is, hasSegment);
1160 }
1161
1166 template<class P = Policy, std::enable_if_t<P::dimIntersection == 0, int> = 0>
1167 static bool intersect_(const Point& a, const Point& b, const Point& c,
1168 const Point& p, const Point& q,
1169 Intersection& is, bool& hasSegmentIs)
1170 {
1171 hasSegmentIs = false;
1172
1173 const auto ab = b - a;
1174 const auto ac = c - a;
1175 const auto qp = p - q;
1176
1177 // compute the triangle normal that defines the triangle plane
1178 const auto normal = crossProduct(ab, ac);
1179
1180 // compute the denominator
1181 // if denom is 0 the segment is parallel and we can return
1182 const auto denom = normal*qp;
1183 const auto abnorm2 = ab.two_norm2();
1184 const auto eps = eps_*abnorm2*qp.two_norm();
1185 using std::abs;
1186 if (abs(denom) < eps)
1187 {
1188 const auto pa = a - p;
1189 const auto denom2 = normal*pa;
1190 if (abs(denom2) > eps_*pa.two_norm()*abnorm2)
1191 return false;
1192
1193 // if we get here, we are in-plane. Check if a
1194 // segment-like intersection with geometry 1 exists.
1195 using SegmentPolicy = typename IntersectionPolicy::SegmentPolicy<ctype, dimworld>;
1196 using Triangle = Dune::AffineGeometry<ctype, 2, dimworld>;
1197 using Segment = Dune::AffineGeometry<ctype, 1, dimworld>;
1198 using SegmentIntersectionAlgorithm = GeometryIntersection<Triangle, Segment, SegmentPolicy>;
1199 using SegmentIntersectionType = typename SegmentIntersectionAlgorithm::Intersection;
1200 SegmentIntersectionType segmentIs;
1201
1202 Triangle triangle(Dune::GeometryTypes::simplex(2), std::array<Point, 3>({a, b, c}));
1203 Segment segment(Dune::GeometryTypes::simplex(1), std::array<Point, 2>({p, q}));
1204 if (SegmentIntersectionAlgorithm::intersection(triangle, segment, segmentIs))
1205 {
1206 hasSegmentIs = true;
1207 return false;
1208 }
1209
1210 // We are in-plane. A point-like
1211 // intersection can only be on the edges
1212 for (const auto& ip : {p, q})
1213 {
1214 if ( intersectsPointSimplex(ip, a, b)
1215 || intersectsPointSimplex(ip, b, c)
1216 || intersectsPointSimplex(ip, c, a) )
1217 {
1218 is = ip;
1219 return true;
1220 }
1221 }
1222
1223 return false;
1224 }
1225
1226 // compute intersection t value of pq with plane of triangle.
1227 // a segment intersects if and only if 0 <= t <= 1.
1228 const auto ap = p - a;
1229 const auto t = (ap*normal)/denom;
1230 if (t < 0.0 - eps_) return false;
1231 if (t > 1.0 + eps_) return false;
1232
1233 // compute the barycentric coordinates and check if the intersection point
1234 // is inside the bounds of the triangle
1235 const auto e = crossProduct(qp, ap);
1236 const auto v = (ac*e)/denom;
1237 if (v < -eps_ || v > 1.0 + eps_) return false;
1238 const auto w = -(ab*e)/denom;
1239 if (w < -eps_ || v + w > 1.0 + eps_) return false;
1240
1241 // Now we are sure there is an intersection points
1242 // Perform delayed division compute the last barycentric coordinate component
1243 const auto u = 1.0 - v - w;
1244
1245 Point ip(0.0);
1246 ip.axpy(u, a);
1247 ip.axpy(v, b);
1248 ip.axpy(w, c);
1249 is = ip;
1250 return true;
1251 }
1252};
1253
1258template <class Geometry1, class Geometry2, class Policy>
1259class GeometryIntersection<Geometry1, Geometry2, Policy, 3, 1, 2>
1260: public GeometryIntersection<Geometry2, Geometry1, Policy, 3, 2, 1>
1261{
1263public:
1271 template<class P = Policy>
1272 static bool intersection(const Geometry1& geo1, const Geometry2& geo2, typename Base::Intersection& intersection)
1273 {
1274 return Base::intersection(geo2, geo1, intersection);
1275 }
1276};
1277
1278} // end namespace Dumux
1279
1280# endif
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.
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
make the local view function available whenever we use the grid geometry
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: geometryintersection.hh:94
constexpr double eps
epsilon for checking direction of scvf normals
Definition: test_tpfafvgeometry_nonconforming.cc:44
Policy structure for point-like intersections.
Definition: geometryintersection.hh:42
static constexpr int dimIntersection
Definition: geometryintersection.hh:47
Point Intersection
Definition: geometryintersection.hh:48
ct ctype
Definition: geometryintersection.hh:43
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:44
static constexpr int dimWorld
Definition: geometryintersection.hh:46
Policy structure for segment-like intersections.
Definition: geometryintersection.hh:54
Segment Intersection
Definition: geometryintersection.hh:61
static constexpr int dimIntersection
Definition: geometryintersection.hh:60
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:56
static constexpr int dimWorld
Definition: geometryintersection.hh:59
std::array< Point, 2 > Segment
Definition: geometryintersection.hh:57
ct ctype
Definition: geometryintersection.hh:55
Policy structure for polygon-like intersections.
Definition: geometryintersection.hh:67
static constexpr int dimWorld
Definition: geometryintersection.hh:72
static constexpr int dimIntersection
Definition: geometryintersection.hh:73
ct ctype
Definition: geometryintersection.hh:68
Polygon Intersection
Definition: geometryintersection.hh:74
Dune::FieldVector< ctype, dw > Point
Definition: geometryintersection.hh:69
std::vector< Point > Polygon
Definition: geometryintersection.hh:70
default policy chooser class
Definition: geometryintersection.hh:80
std::conditional_t< isDim==2, PolygonPolicy< ctype, Geometry1::coorddimension >, SegmentPolicy< ctype, Geometry1::coorddimension > > type
Definition: geometryintersection.hh:89
A class for geometry collision detection and intersection calculation The class can be specialized fo...
Definition: geometryintersection.hh:197
typename Policy::ctype ctype
Definition: geometryintersection.hh:201
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:203
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:206
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:209
typename Policy::Point Point
Definition: geometryintersection.hh:202
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:237
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:234
typename Policy::ctype ctype
Definition: geometryintersection.hh:232
typename Policy::Point Point
Definition: geometryintersection.hh:233
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two segments.
Definition: geometryintersection.hh:246
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:369
typename Policy::ctype ctype
Definition: geometryintersection.hh:349
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:354
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:351
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename P::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:393
typename Policy::Point Point
Definition: geometryintersection.hh:350
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:480
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:500
typename Policy::ctype ctype
Definition: geometryintersection.hh:498
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:503
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding two polygons.
Definition: geometryintersection.hh:524
typename Policy::Point Point
Definition: geometryintersection.hh:499
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:649
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:652
typename Policy::Point Point
Definition: geometryintersection.hh:648
typename Policy::ctype ctype
Definition: geometryintersection.hh:647
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:671
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:781
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:801
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:804
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:827
typename Policy::Point Point
Definition: geometryintersection.hh:800
typename Policy::ctype ctype
Definition: geometryintersection.hh:799
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding polygon and convex polyhedron.
Definition: geometryintersection.hh:976
Intersection IntersectionType
Deprecated alias, will be removed after 3.1.
Definition: geometryintersection.hh:999
typename Policy::Intersection Intersection
Definition: geometryintersection.hh:996
typename Policy::Point Point
Definition: geometryintersection.hh:995
static bool intersection(const Point &a, const Point &b, const Point &c, const Point &p, const Point &q, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1091
typename Policy::ctype ctype
Definition: geometryintersection.hh:994
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &intersection)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1016
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, Intersection &is)
Colliding segment and convex polyhedron.
Definition: geometryintersection.hh:1042
static bool intersection(const Geometry1 &geo1, const Geometry2 &geo2, typename Base::Intersection &intersection)
Colliding segment and convex polygon.
Definition: geometryintersection.hh:1272
An axis-aligned bounding box volume hierarchy for dune grids.