15#ifndef dealii_cgal_utilities_h
16#define dealii_cgal_utilities_h
24#ifdef DEAL_II_WITH_CGAL
29# include <boost/hana.hpp>
31# include <CGAL/Cartesian.h>
32# include <CGAL/Complex_2_in_triangulation_3.h>
33# include <CGAL/Exact_predicates_exact_constructions_kernel.h>
34# include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
35# include <CGAL/Kernel_traits.h>
36# include <CGAL/Mesh_complex_3_in_triangulation_3.h>
37# include <CGAL/Mesh_criteria_3.h>
38# include <CGAL/Mesh_triangulation_3.h>
42# include <CGAL/Polygon_mesh_processing/corefinement.h>
44# include <CGAL/Polygon_mesh_processing/measure.h>
45# include <CGAL/Polygon_mesh_processing/remesh.h>
46# include <CGAL/Polygon_mesh_processing/triangulate_faces.h>
47# include <CGAL/Polyhedral_mesh_domain_with_features_3.h>
48# include <CGAL/Simple_cartesian.h>
49# include <CGAL/Surface_mesh.h>
50# include <CGAL/Triangulation_3.h>
51# include <CGAL/boost/graph/copy_face_graph.h>
52# include <CGAL/boost/graph/selection.h>
53# include <CGAL/convex_hull_3.h>
54# include <CGAL/make_mesh_3.h>
55# include <CGAL/make_surface_mesh.h>
60# include <type_traits>
79# ifdef CGAL_CONCURRENT_MESH_3
140 template <
typename C3t3>
190 template <
typename CGALPo
intType>
209 template <
typename CGALTriangulationType>
212 const unsigned int degree);
229 template <
int dim0,
int dim1,
int spacedim>
232 const typename ::Triangulation<dim0, spacedim>::cell_iterator &
cell0,
233 const typename ::Triangulation<dim1, spacedim>::cell_iterator &
cell1,
234 const unsigned int degree,
237 (ReferenceCells::get_hypercube<dim0>()
240 (ReferenceCells::get_hypercube<dim1>()
258 template <
int dim0,
int dim1,
int spacedim>
261 const typename ::Triangulation<dim0, spacedim>::cell_iterator &
cell0,
262 const typename ::Triangulation<dim1, spacedim>::cell_iterator &
cell1,
263 const unsigned int degree,
295 template <
typename CGALPo
intType>
305 template <
typename C3t3>
314 ExcMessage(
"The surface mesh must be closed."));
316 using K =
typename CGAL::Kernel_traits<CGALPointType>::Kernel;
320 using Tr =
typename CGAL::
321 Mesh_triangulation_3<Mesh_domain, CGAL::Default, ConcurrencyTag>::type;
324 CGAL::Polygon_mesh_processing::triangulate_faces(
surface_mesh);
328 CGAL::parameters::facet_angle =
data.facet_angle,
329 CGAL::parameters::facet_distance =
331 CGAL::parameters::cell_radius_edge_ratio =
332 data.cell_radius_edge_ratio,
333 CGAL::parameters::cell_size =
data.cell_size);
337 CGAL::parameters::no_perturb(),
338 CGAL::parameters::no_exude());
343 template <
typename CGALPo
intType>
353 "output_surface_mesh must be empty upon calling this function"));
356 "The input surface_mesh_1 must be a closed surface mesh."));
359 "The input surface_mesh_2 must be a closed surface mesh."));
361 ExcMessage(
"The first CGAL mesh must be triangulated."));
363 ExcMessage(
"The second CGAL mesh must be triangulated."));
398 ExcMessage(
"The boolean operation was not successfully computed."));
403 template <
typename CGALTriangulationType>
406 const unsigned int degree)
409 Assert(CGALTriangulationType::Point::Ambient_dimension::value == 3,
412 ExcMessage(
"The degree of the Quadrature formula is not positive."));
414 constexpr int spacedim =
415 CGALTriangulationType::Point::Ambient_dimension::value;
416 std::vector<std::array<::Point<spacedim>, spacedim + 1>>
419 std::array<::Point<spacedim>, spacedim + 1>
simplex;
422 for (
unsigned int i = 0; i < (spacedim + 1); ++i)
436 template <
int dim0,
int dim1,
int spacedim>
439 const typename ::Triangulation<dim0, spacedim>::cell_iterator &
cell0,
440 const typename ::Triangulation<dim1, spacedim>::cell_iterator &
cell1,
441 const unsigned int degree,
473 CGAL::Polygon_mesh_processing::triangulate_faces(
surface_1);
474 CGAL::Polygon_mesh_processing::triangulate_faces(
surface_2);
477 ExcMessage(
"The surface must be a triangle mesh."));
488 template <
int dim0,
int dim1,
int spacedim>
491 const typename ::Triangulation<dim0, spacedim>::cell_iterator &
cell0,
492 const typename ::Triangulation<dim1, spacedim>::cell_iterator &
cell1,
493 const unsigned int degree,
507 CGAL::Polygon_mesh_processing::triangulate_faces(
surface_1);
508 CGAL::Polygon_mesh_processing::triangulate_faces(
surface_2);
519 tr.insert(dummy.points().begin(), dummy.points().end());
525 template <
typename CGALPo
intType>
533 using Polyhedron = CGAL::Mesh_polyhedron_3<K>::type;
535 using Tr = CGAL::Mesh_triangulation_3<Mesh_domain>::type;
557 if (
data.edge_size > 0 &&
561 CGAL::parameters::facet_size =
data.facet_size,
562 CGAL::parameters::facet_angle =
data.facet_angle,
563 CGAL::parameters::facet_distance =
565 CGAL::parameters::cell_radius_edge_ratio =
566 data.cell_radius_edge_ratio,
567 CGAL::parameters::cell_size =
data.cell_size);
571 CGAL::parameters::no_perturb(),
572 CGAL::parameters::no_exude());
579 CGAL::parameters::facet_angle =
data.facet_angle,
580 CGAL::parameters::facet_distance =
582 CGAL::parameters::cell_radius_edge_ratio =
583 data.cell_radius_edge_ratio,
584 CGAL::parameters::cell_size =
data.cell_size);
588 CGAL::parameters::no_perturb(),
589 CGAL::parameters::no_exude());
607 template <
int spacedim>
617 if constexpr (spacedim == 2)
620 std::swap(vertices[2], vertices[3]);
632 template <
int dim,
int spacedim>
633 std::vector<Point<spacedim>>
635 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
639 const unsigned int n_vertices = cell->n_vertices();
640 Assert((n_vertices == ReferenceCells::get_hypercube<dim>().n_vertices()) ||
641 (n_vertices == ReferenceCells::get_simplex<dim>().n_vertices()),
645 std::copy_n(mapping.get_vertices(cell).begin(),
664 template <
int dim,
int spacedim>
665 std::vector<Point<spacedim>>
667 const typename ::Triangulation<dim, spacedim>::cell_iterator &cell,
672 const unsigned int n_vertices = cell->face(
face_no)->n_vertices();
674 (n_vertices == ReferenceCells::get_hypercube<dim - 1>().n_vertices()) ||
675 (n_vertices == ReferenceCells::get_simplex<dim - 1>().n_vertices()),
679 std::copy_n(mapping.get_vertices(cell,
face_no).begin(),
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
static ::ExceptionBase & ExcNotImplemented()
#define Assert(cond, exc)
static ::ExceptionBase & ExcMessage(std::string arg1)
#define DEAL_II_ASSERT_UNREACHABLE()
std::vector< index_type > data
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)
CGAL::Sequential_tag ConcurrencyTag
::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 >{})
CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt K
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)
void compute_boolean_operation(const CGAL::Surface_mesh< CGALPointType > &surface_mesh_1, const CGAL::Surface_mesh< CGALPointType > &surface_mesh_2, const BooleanOperation &boolean_operation, CGAL::Surface_mesh< CGALPointType > &output_surface_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 Quadrilateral
::VectorizedArray< Number, width > abs(const ::VectorizedArray< Number, width > &)
const ::parallel::distributed::Triangulation< dim, spacedim > * triangulation