From: Marco Feder Date: Thu, 28 Apr 2022 16:44:52 +0000 (+0000) Subject: Added CGALWrappers::to_cgal_mesh() X-Git-Tag: v9.4.0-rc1~284^2~6 X-Git-Url: https://gitweb.dealii.org/cgi-bin/gitweb.cgi?a=commitdiff_plain;h=684ffe8430e840a32d758d2164c6ab9e7921db3f;p=dealii.git Added CGALWrappers::to_cgal_mesh() --- diff --git a/cmake/config/template-arguments.in b/cmake/config/template-arguments.in index f00dcc1442..4a9cc4aab1 100644 --- a/cmake/config/template-arguments.in +++ b/cmake/config/template-arguments.in @@ -277,3 +277,8 @@ SYM_RANKS := { 2; 4 } OUTPUT_FLAG_TYPES := { DXFlags; UcdFlags; GnuplotFlags; PovrayFlags; EpsFlags; GmvFlags; TecplotFlags; VtkFlags; SvgFlags; Deal_II_IntermediateFlags } + +// CGAL Kernels +CGAL_KERNELS := {CGAL::Simple_cartesian; CGAL::Exact_predicates_exact_constructions_kernel; + CGAL::Exact_predicates_exact_constructions_kernel_with_sqrt; + CGAL::Exact_predicates_inexact_constructions_kernel } diff --git a/doc/news/changes/minor/20220428MarcoFederLucaHeltai b/doc/news/changes/minor/20220428MarcoFederLucaHeltai new file mode 100644 index 0000000000..1b6bdd0242 --- /dev/null +++ b/doc/news/changes/minor/20220428MarcoFederLucaHeltai @@ -0,0 +1,4 @@ +New: Added function CGALWrappers::to_cgal_mesh() to convert a deal.II cell +into a CGAL::Surface_mesh. +
+(Marco Feder, Luca Heltai, 2022/04/28) diff --git a/include/deal.II/cgal/surface_mesh.h b/include/deal.II/cgal/surface_mesh.h new file mode 100644 index 0000000000..c4752ac2b6 --- /dev/null +++ b/include/deal.II/cgal/surface_mesh.h @@ -0,0 +1,76 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2020 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + +#ifndef dealii_cgal_surface_mesh_h +#define dealii_cgal_surface_mesh_h + +#include + +#include + +#include + +#include + +#ifdef DEAL_II_WITH_CGAL +# include + +DEAL_II_NAMESPACE_OPEN + +namespace CGALWrappers +{ + /** + * Build a CGAL::Surface_mesh from a deal.II cell. + * + * The class Surface_mesh implements a halfedge data structure and can be used + * to represent polyhedral surfaces. It is an edge-centered data structure + * capable of maintaining incidence information of vertices, edges, and faces. + * Each edge is represented by two halfedges with opposite orientation. The + * orientation of a face is chosen so that the halfedges around a face are + * oriented counterclockwise. + * + * More information on this class is available at + * https://doc.cgal.org/latest/Surface_mesh/index.html + * + * The function will throw an exception in dimension one. In dimension two it + * generates a surface mesh of the quadrilateral cell or of the triangle cell, + * while in dimension three it will generate the surface mesh of the cell, + * i.e., a polyhedral mesh containing the faces of the input cell. + * + * The generated mesh is useful when performing geometric operations using + * CGAL::Polygon_mesh_processing, i.e., to compute boolean operations on + * cells, splitting, cutting, slicing, etc. + * + * For examples on how to use the resulting CGAL::Surface_mesh see + * https://doc.cgal.org/latest/Polygon_mesh_processing/ + * + * @param[in] cell The input deal.II cell iterator + * @param[in] mapping The mapping used to map the vertices of the cell + * @param[out] mesh The output CGAL::Surface_mesh + */ + template + void + to_cgal_mesh( + const typename dealii::Triangulation::cell_iterator &cell, + const dealii::Mapping & mapping, + CGAL::Surface_mesh & mesh); +} // namespace CGALWrappers + + + +DEAL_II_NAMESPACE_CLOSE + +#endif +#endif diff --git a/include/deal.II/cgal/utilities.h b/include/deal.II/cgal/utilities.h index b0c3020938..930ed2247a 100644 --- a/include/deal.II/cgal/utilities.h +++ b/include/deal.II/cgal/utilities.h @@ -20,8 +20,12 @@ #ifdef DEAL_II_WITH_CGAL # include +# include +# include +# include # include + DEAL_II_NAMESPACE_OPEN /** * Interface to the Computational Geometry Algorithm Library (CGAL). diff --git a/source/CMakeLists.txt b/source/CMakeLists.txt index fed0ee3dc0..09547929b5 100644 --- a/source/CMakeLists.txt +++ b/source/CMakeLists.txt @@ -56,6 +56,7 @@ ADD_SUBDIRECTORY(fe) ADD_SUBDIRECTORY(dofs) ADD_SUBDIRECTORY(lac) ADD_SUBDIRECTORY(base) +ADD_SUBDIRECTORY(cgal) ADD_SUBDIRECTORY(gmsh) ADD_SUBDIRECTORY(grid) ADD_SUBDIRECTORY(hp) diff --git a/source/cgal/CMakeLists.txt b/source/cgal/CMakeLists.txt new file mode 100644 index 0000000000..e69b6671c1 --- /dev/null +++ b/source/cgal/CMakeLists.txt @@ -0,0 +1,35 @@ +## --------------------------------------------------------------------- +## +## Copyright (C) 2012 - 2022 by the deal.II authors +## +## This file is part of the deal.II library. +## +## The deal.II library is free software; you can use it, redistribute +## it, and/or modify it under the terms of the GNU Lesser General +## Public License as published by the Free Software Foundation; either +## version 2.1 of the License, or (at your option) any later version. +## The full text of the license can be found in the file LICENSE.md at +## the top level directory of deal.II. +## +## --------------------------------------------------------------------- + + +INCLUDE_DIRECTORIES(BEFORE ${CMAKE_CURRENT_BINARY_DIR}) + + +SET(_src +surface_mesh.cc +) + + + +SET(_inst +surface_mesh.inst.in +) + +FILE(GLOB _header + ${CMAKE_SOURCE_DIR}/include/deal.II/cgal/*.h + ) + +DEAL_II_ADD_LIBRARY(obj_cgal OBJECT ${_src} ${_header} ${_inst}) +EXPAND_INSTANTIATIONS(obj_cgal "${_inst}") diff --git a/source/cgal/surface_mesh.cc b/source/cgal/surface_mesh.cc new file mode 100644 index 0000000000..8a47384b37 --- /dev/null +++ b/source/cgal/surface_mesh.cc @@ -0,0 +1,121 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2020 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +#include + +#include + +#ifdef DEAL_II_WITH_CGAL + +DEAL_II_NAMESPACE_OPEN + +namespace +{ + template + void + add_facet(const dealiiFace &face, + const Container & deal2cgal, + CGAL_Mesh & mesh, + const bool clockwise_ordering = true) + { + const unsigned nv = face->n_vertices(); + std::vector indices; + + + switch (nv) + { + case 2: + mesh.add_edge(deal2cgal.at(face->vertex_index(0)), + deal2cgal.at(face->vertex_index(1))); + break; + case 3: + indices = {deal2cgal.at(face->vertex_index(0)), + deal2cgal.at(face->vertex_index(1)), + deal2cgal.at(face->vertex_index(2))}; + break; + case 4: + indices = {deal2cgal.at(face->vertex_index(0)), + deal2cgal.at(face->vertex_index(1)), + deal2cgal.at(face->vertex_index(3)), + deal2cgal.at(face->vertex_index(2))}; + break; + default: + Assert(false, ExcInternalError()); + break; + } + auto f = mesh.null_face(); + if (clockwise_ordering) + f = mesh.add_face(indices); + else + { + std::reverse(indices.begin(), indices.end()); + f = mesh.add_face(indices); + } + Assert(f != mesh.null_face(), + ExcInternalError("While trying to build a CGAL facet, " + "CGAL encountered a orientation problem that it " + "was not able to solve.")); + } +} // namespace + + + +# ifndef DOXYGEN +// Template implementations +namespace CGALWrappers +{ + template + void + to_cgal_mesh(const typename Triangulation::cell_iterator &cell, + const Mapping & mapping, + CGAL::Surface_mesh &mesh) + { + Assert(dim > 1, ExcImpossibleInDim(dim)); + using Mesh = CGAL::Surface_mesh; + using Vertex_index = typename Mesh::Vertex_index; + const auto &vertices = mapping.get_vertices(cell); + + std::map deal2cgal; + // Add all vertices to the mesh + // Store CGAL ordering + for (const auto &i : cell->vertex_indices()) + deal2cgal[cell->vertex_index(i)] = + mesh.add_vertex(CGALWrappers::to_cgal(vertices[i])); + + // Add faces + if (dim < 3) + // simplices and quads + add_facet(cell, deal2cgal, mesh); + else + // faces of 3d cells + for (const auto &f : cell->face_indices()) + add_facet(cell->face(f), + deal2cgal, + mesh, + (f % 2 == 0 || cell->n_vertices() != 8)); + } + // explicit instantiations +# include "surface_mesh.inst" + + +} // namespace CGALWrappers +# endif + + +DEAL_II_NAMESPACE_CLOSE + +#endif diff --git a/source/cgal/surface_mesh.inst.in b/source/cgal/surface_mesh.inst.in new file mode 100644 index 0000000000..9910460881 --- /dev/null +++ b/source/cgal/surface_mesh.inst.in @@ -0,0 +1,26 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2022 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + + + +for (dim : DIMENSIONS; spacedim : SPACE_DIMENSIONS; cgal_kernel : CGAL_KERNELS) + { +#if dim <= spacedim + template void to_cgal_mesh( + const typename Triangulation::cell_iterator &cell, + const Mapping & mapping, + CGAL::Surface_mesh & mesh); +#endif + } diff --git a/tests/cgal/cgal_surface_mesh_01.cc b/tests/cgal/cgal_surface_mesh_01.cc new file mode 100644 index 0000000000..267945ff85 --- /dev/null +++ b/tests/cgal/cgal_surface_mesh_01.cc @@ -0,0 +1,72 @@ +// --------------------------------------------------------------------- +// +// Copyright (C) 2020 by the deal.II authors +// +// This file is part of the deal.II library. +// +// The deal.II library is free software; you can use it, redistribute +// it, and/or modify it under the terms of the GNU Lesser General +// Public License as published by the Free Software Foundation; either +// version 2.1 of the License, or (at your option) any later version. +// The full text of the license can be found in the file LICENSE.md at +// the top level directory of deal.II. +// +// --------------------------------------------------------------------- + +// Convert a deal.II cell to a cgal Surface_mesh + +#include + +#include + +#include + +#include +#include + +#include +#include + +#include "../tests.h" + + +using namespace CGALWrappers; +using CGALPoint = CGAL::Point_3>; + +template +void +test() +{ + deallog << "dim= " << dim << ",\t spacedim= " << spacedim << std::endl; + std::vector> d2t = {{}, {2}, {3, 4}, {4, 5, 6, 8}}; + for (const auto nv : d2t[dim]) + { + Triangulation tria; + CGAL::Surface_mesh mesh; + + const auto ref = ReferenceCell::n_vertices_to_type(dim, nv); + const auto mapping = ref.template get_default_mapping(1); + GridGenerator::reference_cell(tria, ref); + + const auto cell = tria.begin_active(); + to_cgal_mesh(cell, *mapping, mesh); + + Assert(mesh.is_valid(), dealii::ExcMessage("The CGAL mesh is not valid")); + deallog << "deal vertices: " << nv << ", cgal vertices " + << mesh.num_vertices() << std::endl; + deallog << "deal faces: " << cell->n_faces() << ", cgal faces " + << mesh.num_faces() << std::endl; + deallog << "Valid mesh: " << std::boolalpha << mesh.is_valid() + << std::endl; + deallog << mesh << std::endl; + } +} + +int +main() +{ + initlog(); + test<2, 2>(); + test<2, 3>(); + test<3, 3>(); +} diff --git a/tests/cgal/cgal_surface_mesh_01.output b/tests/cgal/cgal_surface_mesh_01.output new file mode 100644 index 0000000000..df27e769bd --- /dev/null +++ b/tests/cgal/cgal_surface_mesh_01.output @@ -0,0 +1,113 @@ + +DEAL::dim= 2, spacedim= 2 +DEAL::deal vertices: 3, cgal vertices 3 +DEAL::deal faces: 3, cgal faces 1 +DEAL::Valid mesh: true +DEAL::OFF +3 1 0 +0.00000 0.00000 0.00000 +1.00000 0.00000 0.00000 +0.00000 1.00000 0.00000 +3 0 1 2 + +DEAL::deal vertices: 4, cgal vertices 4 +DEAL::deal faces: 4, cgal faces 1 +DEAL::Valid mesh: true +DEAL::OFF +4 1 0 +0.00000 0.00000 0.00000 +1.00000 0.00000 0.00000 +0.00000 1.00000 0.00000 +1.00000 1.00000 0.00000 +4 0 1 3 2 + +DEAL::dim= 2, spacedim= 3 +DEAL::deal vertices: 3, cgal vertices 3 +DEAL::deal faces: 3, cgal faces 1 +DEAL::Valid mesh: true +DEAL::OFF +3 1 0 +0.00000 0.00000 0.00000 +1.00000 0.00000 0.00000 +0.00000 1.00000 0.00000 +3 0 1 2 + +DEAL::deal vertices: 4, cgal vertices 4 +DEAL::deal faces: 4, cgal faces 1 +DEAL::Valid mesh: true +DEAL::OFF +4 1 0 +0.00000 0.00000 0.00000 +1.00000 0.00000 0.00000 +0.00000 1.00000 0.00000 +1.00000 1.00000 0.00000 +4 0 1 3 2 + +DEAL::dim= 3, spacedim= 3 +DEAL::deal vertices: 4, cgal vertices 4 +DEAL::deal faces: 4, cgal faces 4 +DEAL::Valid mesh: true +DEAL::OFF +4 4 0 +0.00000 0.00000 0.00000 +1.00000 0.00000 0.00000 +0.00000 1.00000 0.00000 +0.00000 0.00000 1.00000 +3 0 1 2 +3 1 0 3 +3 0 2 3 +3 2 1 3 + +DEAL::deal vertices: 5, cgal vertices 5 +DEAL::deal faces: 5, cgal faces 5 +DEAL::Valid mesh: true +DEAL::OFF +5 5 0 +-1.00000 -1.00000 0.00000 +1.00000 -1.00000 0.00000 +-1.00000 1.00000 0.00000 +1.00000 1.00000 0.00000 +0.00000 0.00000 1.00000 +4 0 1 3 2 +3 0 2 4 +3 3 1 4 +3 1 0 4 +3 2 3 4 + +DEAL::deal vertices: 6, cgal vertices 6 +DEAL::deal faces: 5, cgal faces 5 +DEAL::Valid mesh: true +DEAL::OFF +6 5 0 +0.00000 0.00000 0.00000 +1.00000 0.00000 0.00000 +0.00000 1.00000 0.00000 +0.00000 0.00000 1.00000 +1.00000 0.00000 1.00000 +0.00000 1.00000 1.00000 +3 1 0 2 +3 3 4 5 +4 0 1 4 3 +4 1 2 5 4 +4 2 0 3 5 + +DEAL::deal vertices: 8, cgal vertices 8 +DEAL::deal faces: 6, cgal faces 6 +DEAL::Valid mesh: true +DEAL::OFF +8 6 0 +0.00000 0.00000 0.00000 +1.00000 0.00000 0.00000 +0.00000 1.00000 0.00000 +1.00000 1.00000 0.00000 +0.00000 0.00000 1.00000 +1.00000 0.00000 1.00000 +0.00000 1.00000 1.00000 +1.00000 1.00000 1.00000 +4 0 2 6 4 +4 5 7 3 1 +4 0 4 5 1 +4 3 7 6 2 +4 0 1 3 2 +4 6 7 5 4 +