21#ifdef DEAL_II_WITH_CGAL
31# include <CGAL/Boolean_set_operations_2.h>
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>
66 using K = CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt;
67 using K_exact = CGAL::Exact_predicates_exact_constructions_kernel;
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>;
108 template <
class T,
class... Types>
112 return std::get_if<T>(v);
115 template <
class T,
class... Types>
119 return boost::get<T>(v);
130 template <
typename TargetVariant>
131 struct Repackage : boost::static_visitor<TargetVariant>
133 template <
typename T>
135 operator()(
const T &t)
const
137 return TargetVariant(t);
147 template <
typename... Types>
148 std::optional<std::variant<Types...>>
149 convert_boost_to_std(
const boost::optional<boost::variant<Types...>> &x)
158 using std_variant = std::variant<Types...>;
159 return boost::apply_visitor(Repackage<std_variant>(), *x);
169 template <
typename... Types>
170 const std::optional<std::variant<Types...>> &
171 convert_boost_to_std(
const std::optional<std::variant<Types...>> &opt)
183 std::list<CDT::Edge> &border)
185 if (start->info().nesting_level != -1)
189 std::list<Face_handle> queue;
190 queue.push_back(start);
191 while (!queue.empty())
195 if (fh->info().nesting_level == -1)
197 fh->info().nesting_level = index;
198 for (
int i = 0; i < 3; i++)
202 if (n->info().nesting_level == -1)
204 if (ct.is_constrained(e))
219 for (CDT::Face_handle f : cdt.all_face_handles())
221 f->info().nesting_level = -1;
223 std::list<CDT::Edge> border;
225 while (!border.empty())
227 CDT::Edge e = border.front();
230 if (n->info().nesting_level == -1)
232 mark_domains(cdt, n, e.first->info().nesting_level + 1, border);
249 std::vector<CGALPoint2>>>
257 std::array<CGALPoint2, 3> pts0, pts1;
259 std::transform(triangle0.begin(),
262 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
264 std::transform(triangle1.begin(),
267 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
271 return convert_boost_to_std(
272 CGAL::intersection(cgal_triangle0, cgal_triangle1));
276 std::optional<std::variant<CGALPoint2, CGALSegment2>>
284 std::array<CGALPoint2, 3> pts0;
285 std::array<CGALPoint2, 2> pts1;
287 std::transform(triangle.begin(),
290 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
292 std::transform(segment.begin(),
295 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
299 return convert_boost_to_std(
300 CGAL::intersection(cgal_segment, cgal_triangle));
306 std::vector<Polygon_with_holes_2>
313 std::array<CGALPoint2, 4> pts0, pts1;
315 std::transform(rectangle0.begin(),
318 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
320 std::transform(rectangle1.begin(),
323 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
325 const CGALPolygon first_poly{pts0.begin(), pts0.end()};
326 const CGALPolygon second_poly{pts1.begin(), pts1.end()};
328 std::vector<Polygon_with_holes_2> poly_list;
329 CGAL::intersection(first_poly,
331 std::back_inserter(poly_list));
337 std::optional<std::variant<CGALPoint3, CGALSegment3>>
342# if DEAL_II_CGAL_VERSION_GTE(5, 5, 0)
347 std::array<CGALPoint3, 4> pts0;
348 std::array<CGALPoint3, 2> pts1;
350 std::transform(tetrahedron.begin(),
353 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3, 3>);
355 std::transform(segment.begin(),
358 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3, 3>);
360 CGALTetra cgal_tetrahedron{pts0[0], pts0[1], pts0[2], pts0[3]};
362 return convert_boost_to_std(
363 CGAL::intersection(cgal_segment, cgal_tetrahedron));
368 "This function requires a version of CGAL greater or equal than 5.5."));
380 std::vector<CGALPoint3>>>
385# if DEAL_II_CGAL_VERSION_GTE(5, 5, 0)
390 std::array<CGALPoint3, 4> pts0;
391 std::array<CGALPoint3, 3> pts1;
393 std::transform(tetrahedron.begin(),
396 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3, 3>);
398 std::transform(triangle.begin(),
401 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3, 3>);
403 CGALTetra cgal_tetrahedron{pts0[0], pts0[1], pts0[2], pts0[3]};
405 return convert_boost_to_std(
406 CGAL::intersection(cgal_triangle, cgal_tetrahedron));
412 "This function requires a version of CGAL greater or equal than 5.5."));
420 std::vector<std::array<Point<2>, 3>>
428 const auto intersection_test =
431 if (!intersection_test.empty())
433 const auto &poly = intersection_test[0].outer_boundary();
434 const unsigned int size_poly = poly.size();
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>(
445 else if (size_poly >= 4)
448 std::vector<std::array<Point<2>, 3>> collection;
451 cdt.insert_constraint(poly.vertices_begin(),
459 if (f->info().in_domain() &&
460 CGAL::to_double(cdt.triangle(f).area()) > tol)
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))}});
486 std::vector<std::array<Point<2>, 2>>
494 std::array<CGALPoint2, 4> pts;
496 std::transform(quad.begin(),
499 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint2, 2>);
504 CGALWrappers::dealii_point_to_cgal_point<CGALPoint2>(line[0]),
505 CGALWrappers::dealii_point_to_cgal_point<CGALPoint2>(line[1]));
507 cdt.insert_constraint(poly.vertices_begin(), poly.vertices_end(),
true);
508 std::vector<std::array<Point<2>, 2>> vertices;
512 if (f->info().in_domain() &&
513 CGAL::to_double(cdt.triangle(f).area()) > tol &&
514 CGAL::do_intersect(segm, cdt.triangle(f)))
516 const auto intersection =
517 CGAL::intersection(segm, cdt.triangle(f));
521 {{CGALWrappers::cgal_point_to_dealii_point<2>((*s)[0]),
522 CGALWrappers::cgal_point_to_dealii_point<2>((*s)[1])}});
531 std::vector<std::array<Point<3>, 2>>
536# if DEAL_II_CGAL_VERSION_GTE(5, 5, 0)
541 std::array<CGALPoint3_exact, 8> pts;
547 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
550 CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact>(line[0]),
551 CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact>(line[1]));
555 std::vector<std::array<Point<3>, 2>> vertices;
557 cgal_triangulation.insert(pts.begin(), pts.end());
558 for (
const auto &c : cgal_triangulation.finite_cell_handles())
560 const auto &cgal_tetrahedron = cgal_triangulation.tetrahedron(c);
561 if (CGAL::do_intersect(cgal_segment, cgal_tetrahedron))
563 const auto intersection =
564 CGAL::intersection(cgal_segment, cgal_tetrahedron);
568 if (s->squared_length() > tol * tol)
571 {{CGALWrappers::cgal_point_to_dealii_point<3>(
573 CGALWrappers::cgal_point_to_dealii_point<3>(
584 "This function requires a version of CGAL greater or equal than 5.5."));
592 std::vector<std::array<Point<3>, 3>>
597# if DEAL_II_CGAL_VERSION_GTE(5, 5, 0)
602 std::array<CGALPoint3_exact, 8> pts_hex;
603 std::array<CGALPoint3_exact, 4> pts_quad;
609 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
615 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
618 std::vector<std::array<Point<3>, 3>> vertices;
620 triangulation_hexa.insert(pts_hex.begin(), pts_hex.end());
624 triangulation_quad.insert(pts_quad.begin(), pts_quad.end());
626 for (
const auto &c : triangulation_hexa.finite_cell_handles())
628 const auto &tet = triangulation_hexa.tetrahedron(c);
630 for (
const auto &f : triangulation_quad.finite_facets())
632 if (CGAL::do_intersect(tet, triangulation_quad.triangle(f)))
634 const auto intersection =
635 CGAL::intersection(triangulation_quad.triangle(f), tet);
640 if (CGAL::to_double(t->squared_area()) > tol * tol)
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])}});
649 if (
const std::vector<CGALPoint3_exact> *vps =
650 get_if_<std::vector<CGALPoint3_exact>>(&*intersection))
653 tria_inter.insert(vps->begin(), vps->end());
655 for (
auto it = tria_inter.finite_facets_begin();
656 it != tria_inter.finite_facets_end();
659 const auto triangle = tria_inter.triangle(*it);
660 if (CGAL::to_double(triangle.squared_area()) >
663 std::array<Point<3>, 3> verts = {
664 {CGALWrappers::cgal_point_to_dealii_point<3>(
666 CGALWrappers::cgal_point_to_dealii_point<3>(
668 CGALWrappers::cgal_point_to_dealii_point<3>(
671 vertices.push_back(verts);
684 "This function requires a version of CGAL greater or equal than 5.5."));
692 std::vector<std::array<Point<3>, 4>>
700 std::array<CGALPoint3_exact, 8> pts_hex0;
701 std::array<CGALPoint3_exact, 8> pts_hex1;
707 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
713 &CGALWrappers::dealii_point_to_cgal_point<CGALPoint3_exact, 3>);
717 std::vector<std::array<Point<3>, 4>> vertices;
720 tria0.insert(pts_hex0.begin(), pts_hex0.end());
721 tria1.insert(pts_hex1.begin(), pts_hex1.end());
723 for (
const auto &c0 : tria0.finite_cell_handles())
725 const auto &tet0 = tria1.tetrahedron(c0);
726 const auto &tetg0 = CGAL::make_tetrahedron(tet0.vertex(0),
732 for (
const auto &c1 : tria1.finite_cell_handles())
734 const auto &tet1 = tria1.tetrahedron(c1);
735 const auto &tetg1 = CGAL::make_tetrahedron(tet1.vertex(0),
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)
748 triangulation_hexa.insert(sm.points().begin(),
750 for (
const auto &c : triangulation_hexa.finite_cell_handles())
752 const auto &tet = triangulation_hexa.tetrahedron(c);
754 {{CGALWrappers::cgal_point_to_dealii_point<3>(
756 CGALWrappers::cgal_point_to_dealii_point<3>(
758 CGALWrappers::cgal_point_to_dealii_point<3>(
760 CGALWrappers::cgal_point_to_dealii_point<3>(
775 template <
int structdim0,
int structdim1,
int spacedim>
776 std::vector<std::array<Point<spacedim>, structdim1 + 1>>
782 const unsigned int n_vertices0 = vertices0.size();
783 const unsigned int n_vertices1 = vertices1.size();
786 n_vertices0 > 0 || n_vertices1 > 0,
788 "The intersection cannot be computed as at least one of the two cells has no vertices."));
790 if constexpr (structdim0 == 2 && structdim1 == 2 && spacedim == 2)
792 if (n_vertices0 == 4 && n_vertices1 == 4)
799 else if constexpr (structdim0 == 2 && structdim1 == 1 && spacedim == 2)
801 if (n_vertices0 == 4 && n_vertices1 == 2)
808 else if constexpr (structdim0 == 3 && structdim1 == 1 && spacedim == 3)
810 if (n_vertices0 == 8 && n_vertices1 == 2)
817 else if constexpr (structdim0 == 3 && structdim1 == 2 && spacedim == 3)
819 if (n_vertices0 == 8 && n_vertices1 == 4)
826 else if constexpr (structdim0 == 3 && structdim1 == 3 && spacedim == 3)
828 if (n_vertices0 == 8 && n_vertices1 == 8)
845 template <
int structdim0,
int structdim1,
int spacedim>
846 std::vector<std::array<Point<spacedim>, structdim1 + 1>>
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);
867 vertices0, vertices1, tol);
877# include "cgal/intersections.inst"
888template <
int structdim0,
893std::vector<std::array<Point<spacedim>, structdim1 + 1>>
905template <
int structdim0,
int structdim1,
int spacedim>
906std::vector<std::array<Point<spacedim>, structdim1 + 1>>
Abstract base class for mapping classes.
virtual boost::container::small_vector< Point< spacedim >, ReferenceCells::max_n_vertices< dim >() > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_DISABLE_EXTRA_DIAGNOSTICS
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ENABLE_EXTRA_DIAGNOSTICS
#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
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)
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
constexpr const ReferenceCell & get_hypercube()