deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
intersections.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2022 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#include <deal.II/base/config.h>
16
18
19#include <algorithm>
20
21#ifdef DEAL_II_WITH_CGAL
22
25
26# include <deal.II/fe/mapping.h>
27
28# include <deal.II/grid/tria.h>
29
31# include <CGAL/Boolean_set_operations_2.h>
33
35
36# include <CGAL/Cartesian.h>
37# include <CGAL/Circular_kernel_intersections.h>
38# include <CGAL/Constrained_Delaunay_triangulation_2.h>
39# include <CGAL/Delaunay_mesh_face_base_2.h>
40# include <CGAL/Delaunay_mesh_size_criteria_2.h>
41# include <CGAL/Delaunay_mesher_2.h>
42# include <CGAL/Delaunay_triangulation_2.h>
43# include <CGAL/Exact_predicates_exact_constructions_kernel_with_sqrt.h>
44# include <CGAL/Kernel_traits.h>
45# include <CGAL/Polygon_2.h>
46# include <CGAL/Polygon_with_holes_2.h>
47# include <CGAL/Projection_traits_xy_3.h>
48# include <CGAL/Segment_3.h>
49# include <CGAL/Simple_cartesian.h>
50# include <CGAL/Surface_mesh/Surface_mesh.h>
51# include <CGAL/Tetrahedron_3.h>
52# include <CGAL/Triangle_2.h>
53# include <CGAL/Triangle_3.h>
54# include <CGAL/Triangulation_2.h>
55# include <CGAL/Triangulation_3.h>
56# include <CGAL/Triangulation_face_base_with_id_2.h>
57# include <CGAL/Triangulation_face_base_with_info_2.h>
58
59# include <optional>
60# include <variant>
61
63
64namespace CGALWrappers
65{
66 using K = CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt;
67 using K_exact = CGAL::Exact_predicates_exact_constructions_kernel;
68 using CGALPolygon = CGAL::Polygon_2<K>;
69 using Polygon_with_holes_2 = CGAL::Polygon_with_holes_2<K>;
70 using CGALTriangle2 = K::Triangle_2;
71 using CGALTriangle3 = K::Triangle_3;
72 using CGALTriangle3_exact = K_exact::Triangle_3;
73 using CGALPoint2 = K::Point_2;
74 using CGALPoint3 = K::Point_3;
75 using CGALPoint3_exact = K_exact::Point_3;
76 using CGALSegment2 = K::Segment_2;
77 using Surface_mesh = CGAL::Surface_mesh<K_exact::Point_3>;
78 using CGALSegment3 = K::Segment_3;
79 using CGALSegment3_exact = K_exact::Segment_3;
80 using CGALTetra = K::Tetrahedron_3;
81 using CGALTetra_exact = K_exact::Tetrahedron_3;
82 using Triangulation2 = CGAL::Triangulation_2<K>;
83 using Triangulation3 = CGAL::Triangulation_3<K>;
84 using Triangulation3_exact = CGAL::Triangulation_3<K_exact>;
85
86 struct FaceInfo2
87 {
88 FaceInfo2() = default;
90 bool
92 {
93 return nesting_level % 2 == 1;
94 }
95 };
96
97 using Vb = CGAL::Triangulation_vertex_base_2<K>;
98 using Fbb = CGAL::Triangulation_face_base_with_info_2<FaceInfo2, K>;
99 using CFb = CGAL::Constrained_triangulation_face_base_2<K, Fbb>;
100 using Fb = CGAL::Delaunay_mesh_face_base_2<K, CFb>;
101 using Tds = CGAL::Triangulation_data_structure_2<Vb, Fb>;
102 using Itag = CGAL::Exact_predicates_tag;
103 using CDT = CGAL::Constrained_Delaunay_triangulation_2<K, Tds, Itag>;
104 using Criteria = CGAL::Delaunay_mesh_size_criteria_2<CDT>;
105 using Vertex_handle = CDT::Vertex_handle;
106 using Face_handle = CDT::Face_handle;
107
108 template <class T, class... Types>
109 const T *
110 get_if_(const std::variant<Types...> *v)
111 {
112 return std::get_if<T>(v);
113 }
114
115 template <class T, class... Types>
116 const T *
117 get_if_(const boost::variant<Types...> *v)
118 {
119 return boost::get<T>(v);
120 }
121
122 namespace internal
123 {
124 namespace
125 {
130 template <typename TargetVariant>
131 struct Repackage : boost::static_visitor<TargetVariant>
132 {
133 template <typename T>
134 TargetVariant
135 operator()(const T &t) const
136 {
137 return TargetVariant(t);
138 }
139 };
140
147 template <typename... Types>
148 std::optional<std::variant<Types...>>
149 convert_boost_to_std(const boost::optional<boost::variant<Types...>> &x)
150 {
151 if (x)
152 {
153 // The boost::optional object contains an object of type
154 // boost::variant. We need to unpack which type the
155 // variant contains, and re-package that into a
156 // std::variant. This is easily done using a visitor
157 // object.
158 using std_variant = std::variant<Types...>;
159 return boost::apply_visitor(Repackage<std_variant>(), *x);
160 }
161 else
162 {
163 // The boost::optional object was empty. Return an empty
164 // std::optional object.
165 return {};
166 }
167 }
168
169 template <typename... Types>
170 const std::optional<std::variant<Types...>> &
171 convert_boost_to_std(const std::optional<std::variant<Types...>> &opt)
172 {
173 return opt;
174 }
175 } // namespace
176
177
178
179 void
181 Face_handle start,
182 int index,
183 std::list<CDT::Edge> &border)
184 {
185 if (start->info().nesting_level != -1)
186 {
187 return;
188 }
189 std::list<Face_handle> queue;
190 queue.push_back(start);
191 while (!queue.empty())
192 {
193 Face_handle fh = queue.front();
194 queue.pop_front();
195 if (fh->info().nesting_level == -1)
196 {
197 fh->info().nesting_level = index;
198 for (int i = 0; i < 3; i++)
199 {
200 CDT::Edge e(fh, i);
201 Face_handle n = fh->neighbor(i);
202 if (n->info().nesting_level == -1)
203 {
204 if (ct.is_constrained(e))
205 border.push_back(e);
206 else
207 queue.push_back(n);
208 }
209 }
210 }
211 }
212 }
213
214
215
216 void
218 {
219 for (CDT::Face_handle f : cdt.all_face_handles())
220 {
221 f->info().nesting_level = -1;
222 }
223 std::list<CDT::Edge> border;
224 mark_domains(cdt, cdt.infinite_face(), 0, border);
225 while (!border.empty())
226 {
227 CDT::Edge e = border.front();
228 border.pop_front();
229 Face_handle n = e.first->neighbor(e.second);
230 if (n->info().nesting_level == -1)
231 {
232 mark_domains(cdt, n, e.first->info().nesting_level + 1, border);
233 }
234 }
235 }
236
237 // Collection of utilities that compute intersection between simplices
238 // identified by array of points. The return type is the one of
239 // CGAL::intersection(), i.e. a std::optional<std::variant<>>.
240 // Intersection between 2d and 3d objects and 1d/3d objects are available
241 // only with CGAL versions greater or equal than 5.5, hence the
242 // corresponding functions are guarded by #ifdef directives. All the
243 // signatures follow the convection that the first entity has an intrinsic
244 // dimension higher than the second one.
245
246 std::optional<std::variant<CGALPoint2,
249 std::vector<CGALPoint2>>>
251 const ArrayView<const Point<2>> &triangle0,
252 const ArrayView<const Point<2>> &triangle1)
253 {
254 AssertDimension(triangle0.size(), 3);
255 AssertDimension(triangle0.size(), triangle1.size());
256
257 std::array<CGALPoint2, 3> pts0, pts1;
258
259 std::transform(triangle0.begin(),
260 triangle0.end(),
261 pts0.begin(),
262 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
263
264 std::transform(triangle1.begin(),
265 triangle1.end(),
266 pts1.begin(),
267 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
268
269 CGALTriangle2 cgal_triangle0{pts0[0], pts0[1], pts0[2]};
270 CGALTriangle2 cgal_triangle1{pts1[0], pts1[1], pts1[2]};
271 return convert_boost_to_std(
272 CGAL::intersection(cgal_triangle0, cgal_triangle1));
273 }
274
275
276 std::optional<std::variant<CGALPoint2, CGALSegment2>>
278 const ArrayView<const Point<2>> &triangle,
279 const ArrayView<const Point<2>> &segment)
280 {
281 AssertDimension(triangle.size(), 3);
282 AssertDimension(segment.size(), 2);
283
284 std::array<CGALPoint2, 3> pts0;
285 std::array<CGALPoint2, 2> pts1;
286
287 std::transform(triangle.begin(),
288 triangle.end(),
289 pts0.begin(),
290 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
291
292 std::transform(segment.begin(),
293 segment.end(),
294 pts1.begin(),
295 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
296
297 CGALTriangle2 cgal_triangle{pts0[0], pts0[1], pts0[2]};
298 CGALSegment2 cgal_segment{pts1[0], pts1[1]};
299 return convert_boost_to_std(
300 CGAL::intersection(cgal_segment, cgal_triangle));
301 }
302
303
304
305 // rectangle-rectangle
306 std::vector<Polygon_with_holes_2>
308 const ArrayView<const Point<2>> &rectangle1)
309 {
310 AssertDimension(rectangle0.size(), 4);
311 AssertDimension(rectangle0.size(), rectangle1.size());
312
313 std::array<CGALPoint2, 4> pts0, pts1;
314
315 std::transform(rectangle0.begin(),
316 rectangle0.end(),
317 pts0.begin(),
318 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
319
320 std::transform(rectangle1.begin(),
321 rectangle1.end(),
322 pts1.begin(),
323 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
324
325 const CGALPolygon first_poly{pts0.begin(), pts0.end()};
326 const CGALPolygon second_poly{pts1.begin(), pts1.end()};
327
328 std::vector<Polygon_with_holes_2> poly_list;
329 CGAL::intersection(first_poly,
330 second_poly,
331 std::back_inserter(poly_list));
332 return poly_list;
333 }
334
335
336
337 std::optional<std::variant<CGALPoint3, CGALSegment3>>
339 const ArrayView<const Point<3>> &tetrahedron,
340 const ArrayView<const Point<3>> &segment)
341 {
342# if DEAL_II_CGAL_VERSION_GTE(5, 5, 0)
343
344 AssertDimension(tetrahedron.size(), 4);
345 AssertDimension(segment.size(), 2);
346
347 std::array<CGALPoint3, 4> pts0;
348 std::array<CGALPoint3, 2> pts1;
349
350 std::transform(tetrahedron.begin(),
351 tetrahedron.end(),
352 pts0.begin(),
353 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3, 3>);
354
355 std::transform(segment.begin(),
356 segment.end(),
357 pts1.begin(),
358 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3, 3>);
359
360 CGALTetra cgal_tetrahedron{pts0[0], pts0[1], pts0[2], pts0[3]};
361 CGALSegment3 cgal_segment{pts1[0], pts1[1]};
362 return convert_boost_to_std(
363 CGAL::intersection(cgal_segment, cgal_tetrahedron));
364# else
365 Assert(
366 false,
368 "This function requires a version of CGAL greater or equal than 5.5."));
369 (void)tetrahedron;
370 (void)segment;
371 return {};
372# endif
373 }
374
375
376 // tetra, triangle
377 std::optional<std::variant<CGALPoint3,
380 std::vector<CGALPoint3>>>
382 const ArrayView<const Point<3>> &tetrahedron,
383 const ArrayView<const Point<3>> &triangle)
384 {
385# if DEAL_II_CGAL_VERSION_GTE(5, 5, 0)
386
387 AssertDimension(tetrahedron.size(), 4);
388 AssertDimension(triangle.size(), 3);
389
390 std::array<CGALPoint3, 4> pts0;
391 std::array<CGALPoint3, 3> pts1;
392
393 std::transform(tetrahedron.begin(),
394 tetrahedron.end(),
395 pts0.begin(),
396 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3, 3>);
397
398 std::transform(triangle.begin(),
399 triangle.end(),
400 pts1.begin(),
401 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3, 3>);
402
403 CGALTetra cgal_tetrahedron{pts0[0], pts0[1], pts0[2], pts0[3]};
404 CGALTriangle3 cgal_triangle{pts1[0], pts1[1], pts1[2]};
405 return convert_boost_to_std(
406 CGAL::intersection(cgal_triangle, cgal_tetrahedron));
407# else
408
409 Assert(
410 false,
412 "This function requires a version of CGAL greater or equal than 5.5."));
413 (void)tetrahedron;
414 (void)triangle;
415 return {};
416# endif
417 }
418
419 // quad-quad
420 std::vector<std::array<Point<2>, 3>>
422 const ArrayView<const Point<2>> &quad1,
423 const double tol)
424 {
425 AssertDimension(quad0.size(), 4);
426 AssertDimension(quad0.size(), quad1.size());
427
428 const auto intersection_test =
430
431 if (!intersection_test.empty())
432 {
433 const auto &poly = intersection_test[0].outer_boundary();
434 const unsigned int size_poly = poly.size();
435 if (size_poly == 3)
436 {
437 // intersection is a triangle itself, so directly return its
438 // vertices.
439 return {
440 {{CGALWrappers::cgal_point_to_dealii_point<2>(poly.vertex(0)),
441 CGALWrappers::cgal_point_to_dealii_point<2>(poly.vertex(1)),
442 CGALWrappers::cgal_point_to_dealii_point<2>(
443 poly.vertex(2))}}};
444 }
445 else if (size_poly >= 4)
446 {
447 // intersection is a polygon, need to triangulate it.
448 std::vector<std::array<Point<2>, 3>> collection;
449
450 CDT cdt;
451 cdt.insert_constraint(poly.vertices_begin(),
452 poly.vertices_end(),
453 true);
454
456
457 for (Face_handle f : cdt.finite_face_handles())
458 {
459 if (f->info().in_domain() &&
460 CGAL::to_double(cdt.triangle(f).area()) > tol)
461 {
462 collection.push_back(
463 {{CGALWrappers::cgal_point_to_dealii_point<2>(
464 cdt.triangle(f).vertex(0)),
465 CGALWrappers::cgal_point_to_dealii_point<2>(
466 cdt.triangle(f).vertex(1)),
467 CGALWrappers::cgal_point_to_dealii_point<2>(
468 cdt.triangle(f).vertex(2))}});
469 }
470 }
471 return collection;
472 }
473 else
474 {
475 Assert(false, ExcMessage("The polygon is degenerate."));
476 return {};
477 }
478 }
479 else
480 {
481 return {};
482 }
483 }
484
485 // Specialization for quad \cap line
486 std::vector<std::array<Point<2>, 2>>
488 const ArrayView<const Point<2>> &line,
489 const double tol)
490 {
491 AssertDimension(quad.size(), 4);
492 AssertDimension(line.size(), 2);
493
494 std::array<CGALPoint2, 4> pts;
495
496 std::transform(quad.begin(),
497 quad.end(),
498 pts.begin(),
499 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
500
501 CGALPolygon poly(pts.begin(), pts.end());
502
503 CGALSegment2 segm(
504 CGALWrappers::dealii_point_to_cgal_point<CGALPoint2>(line[0]),
505 CGALWrappers::dealii_point_to_cgal_point<CGALPoint2>(line[1]));
506 CDT cdt;
507 cdt.insert_constraint(poly.vertices_begin(), poly.vertices_end(), true);
508 std::vector<std::array<Point<2>, 2>> vertices;
510 for (Face_handle f : cdt.finite_face_handles())
511 {
512 if (f->info().in_domain() &&
513 CGAL::to_double(cdt.triangle(f).area()) > tol &&
514 CGAL::do_intersect(segm, cdt.triangle(f)))
515 {
516 const auto intersection =
517 CGAL::intersection(segm, cdt.triangle(f));
518 if (const CGALSegment2 *s = get_if_<CGALSegment2>(&*intersection))
519 {
520 vertices.push_back(
521 {{CGALWrappers::cgal_point_to_dealii_point<2>((*s)[0]),
522 CGALWrappers::cgal_point_to_dealii_point<2>((*s)[1])}});
523 }
524 }
525 }
526
527 return vertices;
528 }
529
530 // specialization for hex \cap line
531 std::vector<std::array<Point<3>, 2>>
533 const ArrayView<const Point<3>> &line,
534 const double tol)
535 {
536# if DEAL_II_CGAL_VERSION_GTE(5, 5, 0)
537
538 AssertDimension(hexa.size(), 8);
539 AssertDimension(line.size(), 2);
540
541 std::array<CGALPoint3_exact, 8> pts;
542
543 std::transform(
544 hexa.begin(),
545 hexa.end(),
546 pts.begin(),
547 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
548
549 CGALSegment3_exact cgal_segment(
550 CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact>(line[0]),
551 CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact>(line[1]));
552
553 // Subdivide the hex into tetrahedrons, and intersect each one of them
554 // with the line
555 std::vector<std::array<Point<3>, 2>> vertices;
556 Triangulation3_exact cgal_triangulation;
557 cgal_triangulation.insert(pts.begin(), pts.end());
558 for (const auto &c : cgal_triangulation.finite_cell_handles())
559 {
560 const auto &cgal_tetrahedron = cgal_triangulation.tetrahedron(c);
561 if (CGAL::do_intersect(cgal_segment, cgal_tetrahedron))
562 {
563 const auto intersection =
564 CGAL::intersection(cgal_segment, cgal_tetrahedron);
565 if (const CGALSegment3_exact *s =
566 get_if_<CGALSegment3_exact>(&*intersection))
567 {
568 if (s->squared_length() > tol * tol)
569 {
570 vertices.push_back(
571 {{CGALWrappers::cgal_point_to_dealii_point<3>(
572 s->vertex(0)),
573 CGALWrappers::cgal_point_to_dealii_point<3>(
574 s->vertex(1))}});
575 }
576 }
577 }
578 }
579 return vertices;
580# else
581 Assert(
582 false,
584 "This function requires a version of CGAL greater or equal than 5.5."));
585 (void)hexa;
586 (void)line;
587 (void)tol;
588 return {};
589# endif
590 }
591
592 std::vector<std::array<Point<3>, 3>>
594 const ArrayView<const Point<3>> &quad,
595 const double tol)
596 {
597# if DEAL_II_CGAL_VERSION_GTE(5, 5, 0)
598
599 AssertDimension(hexa.size(), 8);
600 AssertDimension(quad.size(), 4);
601
602 std::array<CGALPoint3_exact, 8> pts_hex;
603 std::array<CGALPoint3_exact, 4> pts_quad;
604
605 std::transform(
606 hexa.begin(),
607 hexa.end(),
608 pts_hex.begin(),
609 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
610
611 std::transform(
612 quad.begin(),
613 quad.end(),
614 pts_quad.begin(),
615 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
616
617 // Subdivide hex into tetrahedrons
618 std::vector<std::array<Point<3>, 3>> vertices;
619 Triangulation3_exact triangulation_hexa;
620 triangulation_hexa.insert(pts_hex.begin(), pts_hex.end());
621
622 // Subdivide quad into triangles
623 Triangulation3_exact triangulation_quad;
624 triangulation_quad.insert(pts_quad.begin(), pts_quad.end());
625
626 for (const auto &c : triangulation_hexa.finite_cell_handles())
627 {
628 const auto &tet = triangulation_hexa.tetrahedron(c);
629
630 for (const auto &f : triangulation_quad.finite_facets())
631 {
632 if (CGAL::do_intersect(tet, triangulation_quad.triangle(f)))
633 {
634 const auto intersection =
635 CGAL::intersection(triangulation_quad.triangle(f), tet);
636
637 if (const CGALTriangle3_exact *t =
638 get_if_<CGALTriangle3_exact>(&*intersection))
639 {
640 if (CGAL::to_double(t->squared_area()) > tol * tol)
641 {
642 vertices.push_back(
643 {{cgal_point_to_dealii_point<3>((*t)[0]),
644 cgal_point_to_dealii_point<3>((*t)[1]),
645 cgal_point_to_dealii_point<3>((*t)[2])}});
646 }
647 }
648
649 if (const std::vector<CGALPoint3_exact> *vps =
650 get_if_<std::vector<CGALPoint3_exact>>(&*intersection))
651 {
652 Triangulation3_exact tria_inter;
653 tria_inter.insert(vps->begin(), vps->end());
654
655 for (auto it = tria_inter.finite_facets_begin();
656 it != tria_inter.finite_facets_end();
657 ++it)
658 {
659 const auto triangle = tria_inter.triangle(*it);
660 if (CGAL::to_double(triangle.squared_area()) >
661 tol * tol)
662 {
663 std::array<Point<3>, 3> verts = {
664 {CGALWrappers::cgal_point_to_dealii_point<3>(
665 triangle[0]),
666 CGALWrappers::cgal_point_to_dealii_point<3>(
667 triangle[1]),
668 CGALWrappers::cgal_point_to_dealii_point<3>(
669 triangle[2])}};
670
671 vertices.push_back(verts);
672 }
673 }
674 }
675 }
676 }
677 }
678
679 return vertices;
680# else
681 Assert(
682 false,
684 "This function requires a version of CGAL greater or equal than 5.5."));
685 (void)hexa;
686 (void)quad;
687 (void)tol;
688 return {};
689# endif
690 }
691
692 std::vector<std::array<Point<3>, 4>>
694 const ArrayView<const Point<3>> &hexa1,
695 const double tol)
696 {
697 AssertDimension(hexa0.size(), 8);
698 AssertDimension(hexa0.size(), hexa1.size());
699
700 std::array<CGALPoint3_exact, 8> pts_hex0;
701 std::array<CGALPoint3_exact, 8> pts_hex1;
702
703 std::transform(
704 hexa0.begin(),
705 hexa0.end(),
706 pts_hex0.begin(),
707 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
708
709 std::transform(
710 hexa1.begin(),
711 hexa1.end(),
712 pts_hex1.begin(),
713 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
714
715 Surface_mesh surf0, surf1, sm;
716 // Subdivide hex into tetrahedrons
717 std::vector<std::array<Point<3>, 4>> vertices;
718 Triangulation3_exact tria0, tria1;
719
720 tria0.insert(pts_hex0.begin(), pts_hex0.end());
721 tria1.insert(pts_hex1.begin(), pts_hex1.end());
722
723 for (const auto &c0 : tria0.finite_cell_handles())
724 {
725 const auto &tet0 = tria1.tetrahedron(c0);
726 const auto &tetg0 = CGAL::make_tetrahedron(tet0.vertex(0),
727 tet0.vertex(1),
728 tet0.vertex(2),
729 tet0.vertex(3),
730 surf0);
731 (void)tetg0; // instead of C++ 17s [[maybe unused]]
732 for (const auto &c1 : tria1.finite_cell_handles())
733 {
734 const auto &tet1 = tria1.tetrahedron(c1);
735 const auto &tetg1 = CGAL::make_tetrahedron(tet1.vertex(0),
736 tet1.vertex(1),
737 tet1.vertex(2),
738 tet1.vertex(3),
739 surf1);
740 (void)tetg1; // instead of C++ 17s [[maybe unused]]
741 namespace PMP = CGAL::Polygon_mesh_processing;
742 const bool test_intersection =
743 PMP::corefine_and_compute_intersection(surf0, surf1, sm);
744 if (PMP::volume(sm) > tol && test_intersection)
745 {
746 // Collect tetrahedrons
747 Triangulation3_exact triangulation_hexa;
748 triangulation_hexa.insert(sm.points().begin(),
749 sm.points().end());
750 for (const auto &c : triangulation_hexa.finite_cell_handles())
751 {
752 const auto &tet = triangulation_hexa.tetrahedron(c);
753 vertices.push_back(
754 {{CGALWrappers::cgal_point_to_dealii_point<3>(
755 tet.vertex(0)),
756 CGALWrappers::cgal_point_to_dealii_point<3>(
757 tet.vertex(1)),
758 CGALWrappers::cgal_point_to_dealii_point<3>(
759 tet.vertex(2)),
760 CGALWrappers::cgal_point_to_dealii_point<3>(
761 tet.vertex(3))}});
762 }
763 }
764 surf1.clear();
765 sm.clear();
766 }
767 surf0.clear();
768 }
769 return vertices;
770 }
771
772 } // namespace internal
773
774
775 template <int structdim0, int structdim1, int spacedim>
776 std::vector<std::array<Point<spacedim>, structdim1 + 1>>
778 const ArrayView<const Point<spacedim>> &vertices0,
779 const ArrayView<const Point<spacedim>> &vertices1,
780 const double tol)
781 {
782 const unsigned int n_vertices0 = vertices0.size();
783 const unsigned int n_vertices1 = vertices1.size();
784
785 Assert(
786 n_vertices0 > 0 || n_vertices1 > 0,
788 "The intersection cannot be computed as at least one of the two cells has no vertices."));
789
790 if constexpr (structdim0 == 2 && structdim1 == 2 && spacedim == 2)
791 {
792 if (n_vertices0 == 4 && n_vertices1 == 4)
793 {
795 vertices1,
796 tol);
797 }
798 }
799 else if constexpr (structdim0 == 2 && structdim1 == 1 && spacedim == 2)
800 {
801 if (n_vertices0 == 4 && n_vertices1 == 2)
802 {
804 vertices1,
805 tol);
806 }
807 }
808 else if constexpr (structdim0 == 3 && structdim1 == 1 && spacedim == 3)
809 {
810 if (n_vertices0 == 8 && n_vertices1 == 2)
811 {
813 vertices1,
814 tol);
815 }
816 }
817 else if constexpr (structdim0 == 3 && structdim1 == 2 && spacedim == 3)
818 {
819 if (n_vertices0 == 8 && n_vertices1 == 4)
820 {
822 vertices1,
823 tol);
824 }
825 }
826 else if constexpr (structdim0 == 3 && structdim1 == 3 && spacedim == 3)
827 {
828 if (n_vertices0 == 8 && n_vertices1 == 8)
829 {
831 vertices1,
832 tol);
833 }
834 }
835 else
836 {
838 return {};
839 }
840 (void)tol;
841 return {};
842 }
843
844
845 template <int structdim0, int structdim1, int spacedim>
846 std::vector<std::array<Point<spacedim>, structdim1 + 1>>
850 const Mapping<structdim0, spacedim> &mapping0,
851 const Mapping<structdim1, spacedim> &mapping1,
852 const double tol)
853 {
854 Assert(mapping0.get_vertices(cell0).size() ==
857 Assert(mapping1.get_vertices(cell1).size() ==
860
861 const auto &vertices0 =
862 CGALWrappers::get_vertices_in_cgal_order(cell0, mapping0);
863 const auto &vertices1 =
864 CGALWrappers::get_vertices_in_cgal_order(cell1, mapping1);
865
867 vertices0, vertices1, tol);
868 }
869
870// Explicit instantiations.
871//
872// We don't build the instantiations.inst file if deal.II isn't
873// configured with CGAL, but doxygen doesn't know that and tries to
874// find that file anyway for parsing -- which then of course it fails
875// on. So exclude the following from doxygen consideration.
876# ifndef DOXYGEN
877# include "cgal/intersections.inst"
878# endif
879
880} // namespace CGALWrappers
881
883
884#else
885
887
888template <int structdim0,
889 int structdim1,
890 int spacedim,
891 int n_components0,
892 int n_components1>
893std::vector<std::array<Point<spacedim>, structdim1 + 1>>
895 const std::array<Point<spacedim>, n_components0> &vertices0,
896 const std::array<Point<spacedim>, n_components1> &vertices1,
897 const double tol)
898{
899 (void)vertices0;
900 (void)vertices1;
901 (void)tol;
902 AssertThrow(false, ExcNeedsCGAL());
903}
904
905template <int structdim0, int structdim1, int spacedim>
906std::vector<std::array<Point<spacedim>, structdim1 + 1>>
910 const Mapping<structdim0, spacedim> &mapping0,
911 const Mapping<structdim1, spacedim> &mapping1,
912 const double tol)
913{
914 (void)cell0;
915 (void)cell1;
916 (void)mapping0;
917 (void)mapping1;
918 (void)tol;
919 AssertThrow(false, ExcNeedsCGAL());
920}
921
923
924#endif
Abstract base class for mapping classes.
Definition mapping.h:320
virtual boost::container::small_vector< Point< spacedim >, ReferenceCells::max_n_vertices< dim >() > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
Definition point.h:113
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
Definition config.h:603
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
Definition config.h:647
#define DEAL_II_NOT_IMPLEMENTED()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
#define AssertDimension(dim1, dim2)
static ::ExceptionBase & ExcNeedsCGAL()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
std::vector< std::array< Point< 2 >, 2 > > compute_intersection_quad_line(const ArrayView< const Point< 2 > > &quad, const ArrayView< const Point< 2 > > &line, const double tol)
std::vector< std::array< Point< 3 >, 3 > > compute_intersection_hexa_quad(const ArrayView< const Point< 3 > > &hexa, const ArrayView< const Point< 3 > > &quad, const double tol)
void mark_domains(CDT &ct, Face_handle start, int index, std::list< CDT::Edge > &border)
std::optional< std::variant< CGALPoint3, CGALSegment3 > > compute_intersection_tetra_segment(const ArrayView< const Point< 3 > > &tetrahedron, const ArrayView< const Point< 3 > > &segment)
std::vector< Polygon_with_holes_2 > compute_intersection_rect_rect(const ArrayView< const Point< 2 > > &rectangle0, const ArrayView< const Point< 2 > > &rectangle1)
std::vector< std::array< Point< 3 >, 4 > > compute_intersection_hexa_hexa(const ArrayView< const Point< 3 > > &hexa0, const ArrayView< const Point< 3 > > &hexa1, const double tol)
std::vector< std::array< Point< 3 >, 2 > > compute_intersection_hexa_line(const ArrayView< const Point< 3 > > &hexa, const ArrayView< const Point< 3 > > &line, const double tol)
std::optional< std::variant< CGALPoint3, CGALSegment3, CGALTriangle3, std::vector< CGALPoint3 > > > compute_intersection_tetra_triangle(const ArrayView< const Point< 3 > > &tetrahedron, const ArrayView< const Point< 3 > > &triangle)
std::optional< std::variant< CGALPoint2, CGALSegment2 > > compute_intersection_triangle_segment(const ArrayView< const Point< 2 > > &triangle, const ArrayView< const Point< 2 > > &segment)
std::optional< std::variant< CGALPoint2, CGALSegment2, CGALTriangle2, std::vector< CGALPoint2 > > > compute_intersection_triangle_triangle(const ArrayView< const Point< 2 > > &triangle0, const ArrayView< const Point< 2 > > &triangle1)
std::vector< std::array< Point< 2 >, 3 > > compute_intersection_quad_quad(const ArrayView< const Point< 2 > > &quad0, const ArrayView< const Point< 2 > > &quad1, const double tol)
CGAL::Exact_predicates_tag Itag
K::Triangle_2 CGALTriangle2
CGAL::Triangulation_3< K_exact > Triangulation3_exact
CGAL::Constrained_Delaunay_triangulation_2< K, Tds, Itag > CDT
K::Tetrahedron_3 CGALTetra
CGAL::Exact_predicates_exact_constructions_kernel K_exact
const T * get_if_(const std::variant< Types... > *v)
CGAL::Triangulation_vertex_base_2< K > Vb
CGAL::Triangulation_face_base_with_info_2< FaceInfo2, K > Fbb
CGAL::Polygon_2< K > CGALPolygon
K::Segment_3 CGALSegment3
CGAL::Delaunay_mesh_face_base_2< K, CFb > Fb
CGAL::Triangulation_2< K > Triangulation2
K::Segment_2 CGALSegment2
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
CGAL::Triangulation_3< K > Triangulation3
CGAL::Surface_mesh< K_exact::Point_3 > Surface_mesh
std::vector< std::array< Point< spacedim >, structdim1+1 > > compute_intersection_of_cells(const typename Triangulation< structdim0, spacedim >::cell_iterator &cell0, const typename Triangulation< structdim1, spacedim >::cell_iterator &cell1, const Mapping< structdim0, spacedim > &mapping0, const Mapping< structdim1, spacedim > &mapping1, const double tol=1e-9)
K::Point_2 CGALPoint2
CDT::Face_handle Face_handle
K_exact::Tetrahedron_3 CGALTetra_exact
CGAL::Polygon_with_holes_2< K > Polygon_with_holes_2
K_exact::Triangle_3 CGALTriangle3_exact
K_exact::Segment_3 CGALSegment3_exact
CGAL::Delaunay_mesh_size_criteria_2< CDT > Criteria
CGAL::Constrained_triangulation_face_base_2< K, Fbb > CFb
CGAL::Triangulation_data_structure_2< Vb, Fb > Tds
K::Triangle_3 CGALTriangle3
K_exact::Point_3 CGALPoint3_exact
CDT::Vertex_handle Vertex_handle
K::Point_3 CGALPoint3
constexpr const ReferenceCell & get_hypercube()