15#ifndef dealii_cgal_triangulation_h
16#define dealii_cgal_triangulation_h
28#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/Complex_2_in_triangulation_3.h>
38# include <CGAL/IO/facets_in_complex_2_to_triangle_mesh.h>
39# include <CGAL/Implicit_surface_3.h>
40# include <CGAL/Labeled_mesh_domain_3.h>
41# include <CGAL/Mesh_complex_3_in_triangulation_3.h>
42# include <CGAL/Mesh_criteria_3.h>
43# include <CGAL/Mesh_triangulation_3.h>
44# include <CGAL/Polyhedron_3.h>
45# include <CGAL/Surface_mesh.h>
46# include <CGAL/Surface_mesh_default_triangulation_3.h>
47# include <CGAL/Triangulation_2.h>
48# include <CGAL/Triangulation_3.h>
49# include <CGAL/make_mesh_3.h>
50# include <CGAL/make_surface_mesh.h>
86 template <
int spacedim,
typename CGALTriangulation>
88 add_points_to_cgal_triangulation(
const std::vector<Point<spacedim>> &points,
89 CGALTriangulation &triangulation);
122 template <
typename CGALTriangulation,
int dim,
int spacedim>
124 cgal_triangulation_to_dealii_triangulation(
125 const CGALTriangulation &cgal_triangulation,
126 Triangulation<dim, spacedim> &dealii_triangulation);
151 template <
typename CGALTriangulationType,
152 typename CornerIndexType,
153 typename CurveIndexType>
155 cgal_triangulation_to_dealii_triangulation(
156 const CGAL::Mesh_complex_3_in_triangulation_3<CGALTriangulationType,
160 Triangulation<3> &dealii_triangulation);
183 template <
typename CGAL_MeshType>
185 cgal_surface_mesh_to_dealii_triangulation(
const CGAL_MeshType &cgal_mesh,
186 Triangulation<2, 3> &triangulation);
192 template <
int spacedim,
typename CGALTriangulation>
194 add_points_to_cgal_triangulation(
const std::vector<Point<spacedim>> &points,
195 CGALTriangulation &triangulation)
197 Assert(triangulation.is_valid(),
199 "The triangulation you pass to this function should be a valid "
200 "CGAL triangulation."));
202 std::vector<typename CGALTriangulation::Point> cgal_points(points.size());
203 std::transform(points.begin(),
207 return CGALWrappers::dealii_point_to_cgal_point<
208 typename CGALTriangulation::Point>(p);
211 triangulation.insert(cgal_points.begin(), cgal_points.end());
212 Assert(triangulation.is_valid(),
214 "The Triangulation is no longer valid after inserting the points. "
220 template <
typename CGALTriangulationType,
221 typename CornerIndexType,
222 typename CurveIndexType>
224 cgal_triangulation_to_dealii_triangulation(
225 const CGAL::Mesh_complex_3_in_triangulation_3<CGALTriangulationType,
229 Triangulation<3> &dealii_triangulation)
231 using C3T3 = CGAL::Mesh_complex_3_in_triangulation_3<CGALTriangulationType,
236 std::vector<Point<3>> dealii_vertices;
237 std::map<typename C3T3::Vertex_handle, unsigned int>
238 cgal_to_dealii_vertex_map;
240 std::size_t inum = 0;
241 for (
auto vit = cgal_triangulation.triangulation().finite_vertices_begin();
242 vit != cgal_triangulation.triangulation().finite_vertices_end();
245 if (vit->in_dimension() <= -1)
247 cgal_to_dealii_vertex_map[vit] = inum++;
248 dealii_vertices.emplace_back(
249 CGALWrappers::cgal_point_to_dealii_point<3>(vit->point()));
253 std::vector<CellData<3>> cells;
254 for (
auto cgal_cell = cgal_triangulation.cells_in_complex_begin();
255 cgal_cell != cgal_triangulation.cells_in_complex_end();
259 for (
unsigned int i = 0; i < 4; ++i)
260 cell.vertices[i] = cgal_to_dealii_vertex_map[cgal_cell->vertex(i)];
261 cell.manifold_id = cgal_triangulation.subdomain_index(cgal_cell);
262 cells.push_back(cell);
266 SubCellData subcell_data;
267 if constexpr (std::is_integral_v<typename C3T3::Surface_patch_index>)
269 for (
auto face = cgal_triangulation.facets_in_complex_begin();
270 face != cgal_triangulation.facets_in_complex_end();
273 const auto &[cgal_cell, cgal_vertex_face_index] = *face;
279 for (
int i = 0; i < 4; ++i)
280 if (i != cgal_vertex_face_index)
281 dealii_face.vertices[j++] =
282 cgal_to_dealii_vertex_map[cgal_cell->vertex(i)];
283 dealii_face.manifold_id =
284 cgal_triangulation.surface_patch_index(cgal_cell,
285 cgal_vertex_face_index);
290 if constexpr (std::is_integral_v<typename C3T3::Curve_index>)
292 for (
auto edge = cgal_triangulation.edges_in_complex_begin();
293 edge != cgal_triangulation.edges_in_complex_end();
296 const auto &[cgal_cell, v1, v2] = *edge;
298 dealii_edge.vertices[0] =
299 cgal_to_dealii_vertex_map[cgal_cell->vertex(v1)];
300 dealii_edge.vertices[1] =
301 cgal_to_dealii_vertex_map[cgal_cell->vertex(v2)];
302 dealii_edge.manifold_id =
303 cgal_triangulation.curve_index(cgal_cell->vertex(v1),
304 cgal_cell->vertex(v2));
316 template <
typename CGALTriangulation,
int dim,
int spacedim>
318 cgal_triangulation_to_dealii_triangulation(
319 const CGALTriangulation &cgal_triangulation,
320 Triangulation<dim, spacedim> &dealii_triangulation)
323 ExcMessage(
"The dimension of the input CGAL triangulation (" +
324 std::to_string(cgal_triangulation.dimension()) +
325 ") does not match the dimension of the output "
326 "deal.II triangulation (" +
327 std::to_string(dim) +
")."));
330 ExcMessage(
"The output triangulation object needs to be empty."));
333 std::vector<Point<spacedim>> vertices(
334 cgal_triangulation.number_of_vertices());
335 std::vector<CellData<dim>> cells;
336 SubCellData subcell_data;
339 std::map<typename CGALTriangulation::Vertex_handle, unsigned int>
343 for (
const auto &v : cgal_triangulation.finite_vertex_handles())
346 CGALWrappers::cgal_point_to_dealii_point<spacedim>(v->point());
351 const auto has_faces = boost::hana::is_valid(
352 [](
auto &&obj) ->
decltype(obj.finite_face_handles()) {});
354 const auto has_cells = boost::hana::is_valid(
355 [](
auto &&obj) ->
decltype(obj.finite_cell_handles()) {});
358 if constexpr (
decltype(has_faces(cgal_triangulation)){})
361 if (cgal_triangulation.dimension() == 2)
362 for (
const auto &f : cgal_triangulation.finite_face_handles())
365 for (
unsigned int i = 0;
368 cell.vertices[i] = vertex_map[f->vertex(i)];
369 cells.push_back(cell);
371 else if (cgal_triangulation.dimension() == 1)
373 for (
const auto &e : cgal_triangulation.finite_edges())
377 const auto &f =
e.first;
378 const auto &i =
e.second;
386 for (
unsigned int j = 0;
390 cell.vertices[
id++] = vertex_map[f->vertex(j)];
391 cells.push_back(cell);
398 else if constexpr (
decltype(has_cells(cgal_triangulation)){})
401 if (cgal_triangulation.dimension() == 3)
402 for (
const auto &c : cgal_triangulation.finite_cell_handles())
405 for (
unsigned int i = 0;
408 cell.vertices[i] = vertex_map[c->vertex(i)];
409 cells.push_back(cell);
411 else if (cgal_triangulation.dimension() == 2)
413 for (
const auto &facet : cgal_triangulation.finite_facets())
417 const auto &c = facet.first;
418 const auto &i = facet.second;
426 for (
unsigned int j = 0;
430 cell.vertices[
id++] = vertex_map[c->vertex(j)];
431 cells.push_back(cell);
433 else if (cgal_triangulation.dimension() == 1)
435 for (
const auto &edge : cgal_triangulation.finite_edges())
438 const auto &[c, i, j] = edge;
440 cell.vertices[0] = vertex_map[c->vertex(i)];
441 cell.vertices[1] = vertex_map[c->vertex(j)];
442 cells.push_back(cell);
454 template <
typename CGAL_MeshType>
456 cgal_surface_mesh_to_dealii_triangulation(
const CGAL_MeshType &cgal_mesh,
457 Triangulation<2, 3> &triangulation)
461 "Triangulation must be empty upon calling this function."));
463 const auto is_surface_mesh =
464 boost::hana::is_valid([](
auto &&obj) ->
decltype(obj.faces()) {});
466 const auto is_polyhedral =
467 boost::hana::is_valid([](
auto &&obj) ->
decltype(obj.facets_begin()) {});
470 std::vector<::Point<3>> vertices;
471 std::vector<CellData<2>> cells;
472 SubCellData subcell_data;
475 if constexpr (
decltype(is_surface_mesh(cgal_mesh)){})
479 vertices.reserve(cgal_mesh.num_vertices());
480 std::map<typename CGAL_MeshType::Vertex_index, unsigned int> vertex_map;
483 for (
const auto &v : cgal_mesh.vertices())
485 vertices.emplace_back(CGALWrappers::cgal_point_to_dealii_point<3>(
486 cgal_mesh.point(v)));
492 for (
const auto &face : cgal_mesh.faces())
494 const auto face_vertices =
495 CGAL::vertices_around_face(cgal_mesh.halfedge(face), cgal_mesh);
497 AssertThrow(face_vertices.size() == 3 || face_vertices.size() == 4,
498 ExcMessage(
"Only triangle or quadrilateral surface "
499 "meshes are supported in deal.II"));
501 CellData<2> c(face_vertices.size());
502 auto it_vertex = c.vertices.begin();
503 for (
const auto &v : face_vertices)
504 *(it_vertex++) = vertex_map[v];
507 if (face_vertices.size() == 4)
508 std::swap(c.vertices[3], c.vertices[2]);
510 cells.emplace_back(c);
513 else if constexpr (
decltype(is_polyhedral(cgal_mesh)){})
517 vertices.reserve(cgal_mesh.size_of_vertices());
518 std::map<
decltype(cgal_mesh.vertices_begin()),
unsigned int> vertex_map;
521 for (
auto it = cgal_mesh.vertices_begin();
522 it != cgal_mesh.vertices_end();
525 vertices.emplace_back(
526 CGALWrappers::cgal_point_to_dealii_point<3>(it->point()));
527 vertex_map[it] = i++;
532 for (
auto face = cgal_mesh.facets_begin();
533 face != cgal_mesh.facets_end();
536 auto j = face->facet_begin();
539 ExcMessage(
"Only triangle or quadrilateral surface "
540 "meshes are supported in deal.II. You "
541 "tried to read a mesh where a face has " +
543 " vertices per face."));
546 auto it = c.vertices.begin();
549 *(it++) = vertex_map[j->vertex()];
554 std::swap(c.vertices[3], c.vertices[2]);
556 cells.emplace_back(c);
563 "Unsupported CGAL surface triangulation type."));
unsigned int n_vertices() const
virtual void create_triangulation(const std::vector< Point< spacedim > > &vertices, const std::vector< CellData< dim > > &cells, const SubCellData &subcelldata)
unsigned int n_cells() const
#define DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_ASSERT_UNREACHABLE()
constexpr unsigned int GeometryInfo< dim >::vertices_per_face
#define Assert(cond, exc)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
#define AssertThrow(cond, exc)
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
constexpr ReferenceCell Triangle
constexpr ReferenceCell Tetrahedron
constexpr ReferenceCell Line
std::vector< CellData< 2 > > boundary_quads
std::vector< CellData< 1 > > boundary_lines