15#ifndef dealii_cgal_utilities_h
16#define dealii_cgal_utilities_h
24#ifdef DEAL_II_WITH_CGAL
31# include <boost/hana.hpp>
33# include <CGAL/version.h>
34# if CGAL_VERSION_MAJOR >= 6
35# include <CGAL/Installation/internal/disable_deprecation_warnings_and_errors.h>
37# include <CGAL/Cartesian.h>
38# include <CGAL/Complex_2_in_triangulation_3.h>
39# include <CGAL/Exact_predicates_exact_constructions_kernel.h>
40# include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
41# include <CGAL/Kernel_traits.h>
42# include <CGAL/Mesh_complex_3_in_triangulation_3.h>
43# include <CGAL/Mesh_criteria_3.h>
44# include <CGAL/Mesh_triangulation_3.h>
48# include <CGAL/Polygon_mesh_processing/corefinement.h>
50# include <CGAL/Polygon_mesh_processing/measure.h>
51# include <CGAL/Polygon_mesh_processing/remesh.h>
52# include <CGAL/Polygon_mesh_processing/triangulate_faces.h>
53# include <CGAL/Polyhedral_mesh_domain_with_features_3.h>
54# include <CGAL/Simple_cartesian.h>
55# include <CGAL/Surface_mesh.h>
56# include <CGAL/Triangulation_3.h>
57# include <CGAL/boost/graph/copy_face_graph.h>
58# include <CGAL/boost/graph/selection.h>
59# include <CGAL/convex_hull_3.h>
60# include <CGAL/make_mesh_3.h>
61# include <CGAL/make_surface_mesh.h>
65# include <type_traits>
84# ifdef CGAL_CONCURRENT_MESH_3
145 template <
typename C3t3>
148 CGAL::Surface_mesh<typename C3t3::Point::Point> &surface_mesh,
195 template <
typename CGALPo
intType>
198 const CGAL::Surface_mesh<CGALPointType> &surface_mesh_1,
199 const CGAL::Surface_mesh<CGALPointType> &surface_mesh_2,
201 CGAL::Surface_mesh<CGALPointType> &output_surface_mesh);
214 template <
typename CGALTriangulationType>
217 const unsigned int degree);
234 template <
int dim0,
int dim1,
int spacedim>
237 const typename ::Triangulation<dim0, spacedim>::cell_iterator &cell0,
238 const typename ::Triangulation<dim1, spacedim>::cell_iterator &cell1,
239 const unsigned int degree,
263 template <
int dim0,
int dim1,
int spacedim>
266 const typename ::Triangulation<dim0, spacedim>::cell_iterator &cell0,
267 const typename ::Triangulation<dim1, spacedim>::cell_iterator &cell1,
268 const unsigned int degree,
300 template <
typename CGALPo
intType>
310 template <
typename C3t3>
313 CGAL::Surface_mesh<typename C3t3::Point::Point> &surface_mesh,
317 using CGALPointType =
typename C3t3::Point::Point;
318 Assert(CGAL::is_closed(surface_mesh),
319 ExcMessage(
"The surface mesh must be closed."));
321 using K =
typename CGAL::Kernel_traits<CGALPointType>::Kernel;
322 using Mesh_domain = CGAL::Polyhedral_mesh_domain_with_features_3<
324 CGAL::Surface_mesh<CGALPointType>>;
325 using Tr =
typename CGAL::
326 Mesh_triangulation_3<Mesh_domain, CGAL::Default, ConcurrencyTag>::type;
327 using Mesh_criteria = CGAL::Mesh_criteria_3<Tr>;
329 CGAL::Polygon_mesh_processing::triangulate_faces(surface_mesh);
330 Mesh_domain domain(surface_mesh);
331 domain.detect_features();
332 Mesh_criteria criteria(CGAL::parameters::facet_size = data.facet_size,
333 CGAL::parameters::facet_angle = data.facet_angle,
334 CGAL::parameters::facet_distance =
336 CGAL::parameters::cell_radius_edge_ratio =
337 data.cell_radius_edge_ratio,
338 CGAL::parameters::cell_size = data.cell_size);
340 triangulation = CGAL::make_mesh_3<C3t3>(domain,
342 CGAL::parameters::no_perturb(),
343 CGAL::parameters::no_exude());
348 template <
typename CGALPo
intType>
351 const CGAL::Surface_mesh<CGALPointType> &surface_mesh_1,
352 const CGAL::Surface_mesh<CGALPointType> &surface_mesh_2,
354 CGAL::Surface_mesh<CGALPointType> &output_surface_mesh)
356 Assert(output_surface_mesh.is_empty(),
358 "output_surface_mesh must be empty upon calling this function"));
359 Assert(CGAL::is_closed(surface_mesh_1),
361 "The input surface_mesh_1 must be a closed surface mesh."));
362 Assert(CGAL::is_closed(surface_mesh_2),
364 "The input surface_mesh_2 must be a closed surface mesh."));
365 Assert(CGAL::is_triangle_mesh(surface_mesh_1),
366 ExcMessage(
"The first CGAL mesh must be triangulated."));
367 Assert(CGAL::is_triangle_mesh(surface_mesh_2),
368 ExcMessage(
"The second CGAL mesh must be triangulated."));
371 auto surf_1 = surface_mesh_1;
372 auto surf_2 = surface_mesh_2;
373 namespace PMP = CGAL::Polygon_mesh_processing;
374 switch (boolean_operation)
377 res = PMP::corefine_and_compute_union(surf_1,
379 output_surface_mesh);
382 res = PMP::corefine_and_compute_intersection(surf_1,
384 output_surface_mesh);
387 res = PMP::corefine_and_compute_difference(surf_1,
389 output_surface_mesh);
392 PMP::corefine(surf_1,
394 output_surface_mesh = std::move(surf_1);
398 output_surface_mesh.clear();
403 ExcMessage(
"The boolean operation was not successfully computed."));
408 template <
typename CGALTriangulationType>
409 ::Quadrature<CGALTriangulationType::Point::Ambient_dimension::value>
411 const unsigned int degree)
414 Assert(CGALTriangulationType::Point::Ambient_dimension::value == 3,
417 ExcMessage(
"The degree of the Quadrature formula is not positive."));
419 constexpr int spacedim =
420 CGALTriangulationType::Point::Ambient_dimension::value;
421 std::vector<std::array<::Point<spacedim>, spacedim + 1>>
424 std::array<::Point<spacedim>, spacedim + 1>
simplex;
425 for (
const auto &f : tria.finite_cell_handles())
427 for (
unsigned int i = 0; i < (spacedim + 1); ++i)
430 cgal_point_to_dealii_point<spacedim>(f->vertex(i)->point());
433 vec_of_simplices.push_back(simplex);
436 return QGaussSimplex<spacedim>(degree).mapped_quadrature(vec_of_simplices);
441 template <
int dim0,
int dim1,
int spacedim>
442 ::Quadrature<spacedim>
444 const typename ::Triangulation<dim0, spacedim>::cell_iterator &cell0,
445 const typename ::Triangulation<dim1, spacedim>::cell_iterator &cell1,
446 const unsigned int degree,
448 const Mapping<dim0, spacedim> &mapping0,
449 const Mapping<dim1, spacedim> &mapping1)
451 Assert(dim0 == 3 && dim1 == 3 && spacedim == 3,
467 cell0, cell1, degree, mapping0, mapping1);
471 using K = CGAL::Exact_predicates_inexact_constructions_kernel;
472 using CGALPoint = CGAL::Point_3<K>;
473 using CGALTriangulation = CGAL::Triangulation_3<K>;
474 CGAL::Surface_mesh<CGALPoint> surface_1, surface_2, out_surface;
478 CGAL::Polygon_mesh_processing::triangulate_faces(surface_1);
479 CGAL::Polygon_mesh_processing::triangulate_faces(surface_2);
480 Assert(CGAL::is_triangle_mesh(surface_1) &&
481 CGAL::is_triangle_mesh(surface_2),
482 ExcMessage(
"The surface must be a triangle mesh."));
485 CGALTriangulation tria;
486 tria.insert(out_surface.points().begin(), out_surface.points().end());
493 template <
int dim0,
int dim1,
int spacedim>
494 ::Quadrature<spacedim>
496 const typename ::Triangulation<dim0, spacedim>::cell_iterator &cell0,
497 const typename ::Triangulation<dim1, spacedim>::cell_iterator &cell1,
498 const unsigned int degree,
499 const Mapping<dim0, spacedim> &mapping0,
500 const Mapping<dim1, spacedim> &mapping1)
502 Assert(dim0 == 3 && dim1 == 3 && spacedim == 3,
504 using K = CGAL::Exact_predicates_inexact_constructions_kernel;
505 using CGALPoint = CGAL::Point_3<K>;
506 using CGALTriangulation = CGAL::Triangulation_3<K>;
508 CGAL::Surface_mesh<CGALPoint> surface_1, surface_2, out_surface;
512 CGAL::Polygon_mesh_processing::triangulate_faces(surface_1);
513 CGAL::Polygon_mesh_processing::triangulate_faces(surface_2);
519 CGAL::Surface_mesh<CGALPoint> dummy;
520 CGALTriangulation tr;
521 CGAL::convex_hull_3(out_surface.points().begin(),
522 out_surface.points().end(),
524 tr.insert(dummy.points().begin(), dummy.points().end());
530 template <
typename CGALPo
intType>
535 using K = CGAL::Exact_predicates_inexact_constructions_kernel;
536 using Mesh_domain = CGAL::Polyhedral_mesh_domain_with_features_3<K>;
538 using Polyhedron = CGAL::Mesh_polyhedron_3<K>::type;
540 using Tr = CGAL::Mesh_triangulation_3<Mesh_domain>::type;
542 CGAL::Mesh_complex_3_in_triangulation_3<Tr,
543 Mesh_domain::Corner_index,
544 Mesh_domain::Curve_index>;
546 using Mesh_criteria = CGAL::Mesh_criteria_3<Tr>;
549 CGAL::copy_face_graph(cgal_mesh, poly);
552 std::vector<Polyhedron *> poly_ptrs_vector(1, &poly);
556 Mesh_domain domain(poly_ptrs_vector.begin(), poly_ptrs_vector.end());
558 domain.detect_features();
561 const double default_value_edge_size = std::numeric_limits<double>::max();
562 if (data.edge_size > 0 &&
563 std::abs(data.edge_size - default_value_edge_size) > 1e-12)
565 Mesh_criteria criteria(CGAL::parameters::edge_size = data.edge_size,
566 CGAL::parameters::facet_size = data.facet_size,
567 CGAL::parameters::facet_angle = data.facet_angle,
568 CGAL::parameters::facet_distance =
570 CGAL::parameters::cell_radius_edge_ratio =
571 data.cell_radius_edge_ratio,
572 CGAL::parameters::cell_size = data.cell_size);
574 C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain,
576 CGAL::parameters::no_perturb(),
577 CGAL::parameters::no_exude());
579 CGAL::facets_in_complex_3_to_triangle_mesh(c3t3, cgal_mesh);
581 else if (
std::abs(data.edge_size - default_value_edge_size) <= 1e-12)
583 Mesh_criteria criteria(CGAL::parameters::facet_size = data.facet_size,
584 CGAL::parameters::facet_angle = data.facet_angle,
585 CGAL::parameters::facet_distance =
587 CGAL::parameters::cell_radius_edge_ratio =
588 data.cell_radius_edge_ratio,
589 CGAL::parameters::cell_size = data.cell_size);
591 C3t3 c3t3 = CGAL::make_mesh_3<C3t3>(domain,
593 CGAL::parameters::no_perturb(),
594 CGAL::parameters::no_exude());
596 CGAL::facets_in_complex_3_to_triangle_mesh(c3t3, cgal_mesh);
612 template <
int spacedim>
614 resort_dealii_vertices_to_cgal_order(
const unsigned int structdim,
615 std::vector<Point<spacedim>> &vertices)
622 if constexpr (spacedim == 2)
625 std::swap(vertices[2], vertices[3]);
637 template <
int dim,
int spacedim>
638 std::vector<Point<spacedim>>
639 get_vertices_in_cgal_order(
640 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
641 const Mapping<dim, spacedim> &
mapping)
644 const unsigned int n_vertices = cell->n_vertices();
649 std::vector<Point<spacedim>> ordered_vertices(n_vertices);
652 ordered_vertices.begin());
654 resort_dealii_vertices_to_cgal_order(dim, ordered_vertices);
656 return ordered_vertices;
669 template <
int dim,
int spacedim>
670 std::vector<Point<spacedim>>
671 get_vertices_in_cgal_order(
672 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
673 const unsigned int face_no,
674 const Mapping<dim, spacedim> &
mapping)
677 const unsigned int n_vertices = cell->face(face_no)->n_vertices();
683 std::vector<Point<spacedim>> ordered_vertices(n_vertices);
686 ordered_vertices.begin());
688 resort_dealii_vertices_to_cgal_order(dim - 1, ordered_vertices);
690 return ordered_vertices;
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
static ReferenceCell n_vertices_to_type(const int dim, const unsigned int n_vertices)
#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_ASSERT_UNREACHABLE()
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
const Mapping< dim, spacedim > & get_default_linear_mapping(const Triangulation< dim, spacedim > &triangulation)
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
void remesh_surface(CGAL::Surface_mesh< CGALPointType > &surface_mesh, const AdditionalData< 3 > &data=AdditionalData< 3 >{})
::Quadrature< CGALTriangulationType::Point::Ambient_dimension::value > compute_quadrature(const CGALTriangulationType &tria, const unsigned int degree)
std::vector< CGAL::Polygon_with_holes_2< KernelType > > compute_boolean_operation(const CGAL::Polygon_with_holes_2< KernelType > &polygon_1, const CGAL::Polygon_with_holes_2< KernelType > &polygon_2, const BooleanOperation &boolean_operation)
CGAL::Sequential_tag ConcurrencyTag
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
::Quadrature< spacedim > compute_quadrature_on_boolean_operation(const typename ::Triangulation< dim0, spacedim >::cell_iterator &cell0, const typename ::Triangulation< dim1, spacedim >::cell_iterator &cell1, const unsigned int degree, const BooleanOperation &bool_op, const Mapping< dim0, spacedim > &mapping0=(ReferenceCells::get_hypercube< dim0 >() .template get_default_linear_mapping< dim0, spacedim >()), const Mapping< dim1, spacedim > &mapping1=(ReferenceCells::get_hypercube< dim1 >() .template get_default_linear_mapping< dim1, spacedim >()))
void cgal_surface_mesh_to_cgal_triangulation(CGAL::Surface_mesh< typename C3t3::Point::Point > &surface_mesh, C3t3 &triangulation, const AdditionalData< 3 > &data=AdditionalData< 3 >{})
void dealii_cell_to_cgal_surface_mesh(const typename ::Triangulation< dim, spacedim >::cell_iterator &cell, const ::Mapping< dim, spacedim > &mapping, CGAL::Surface_mesh< CGALPointType > &mesh)
::Quadrature< spacedim > compute_quadrature_on_intersection(const typename ::Triangulation< dim0, spacedim >::cell_iterator &cell0, const typename ::Triangulation< dim1, spacedim >::cell_iterator &cell1, const unsigned int degree, const Mapping< dim0, spacedim > &mapping0, const Mapping< dim1, spacedim > &mapping1)
void simplex(Triangulation< dim, dim > &tria, const std::vector< Point< dim > > &vertices)
constexpr const ReferenceCell & get_hypercube()
constexpr ReferenceCell Quadrilateral
constexpr const ReferenceCell & get_simplex()
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)