deal.II version 9.7.0
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
surface_mesh.cc
Go to the documentation of this file.
1// ------------------------------------------------------------------------
2//
3// SPDX-License-Identifier: LGPL-2.1-or-later
4// Copyright (C) 2022 - 2025 by the deal.II authors
5//
6// This file is part of the deal.II library.
7//
8// Part of the source code is dual licensed under Apache-2.0 WITH
9// LLVM-exception OR LGPL-2.1-or-later. Detailed license information
10// governing the source code and code contributions can be found in
11// LICENSE.md and CONTRIBUTING.md at the top level directory of deal.II.
12//
13// ------------------------------------------------------------------------
14
15#include <deal.II/base/config.h>
16
18
19#ifdef DEAL_II_WITH_CGAL
20
21# include <CGAL/Exact_predicates_exact_constructions_kernel.h>
22# include <CGAL/Simple_cartesian.h>
23# include <CGAL/Surface_mesh/Surface_mesh.h>
24
25
27
28namespace
29{
30 template <typename dealiiFace, typename CGAL_Mesh>
31 void
32 add_facet(
33 const dealiiFace &face,
34 const std::map<unsigned int, typename CGAL_Mesh::Vertex_index> &deal2cgal,
35 CGAL_Mesh &mesh,
36 const bool clockwise_ordering = true)
37 {
38 const auto reference_cell_type = face->reference_cell();
39 std::vector<typename CGAL_Mesh::Vertex_index> indices;
40
41 if (reference_cell_type == ReferenceCells::Line)
42 {
43 mesh.add_edge(deal2cgal.at(face->vertex_index(0)),
44 deal2cgal.at(face->vertex_index(1)));
45
46 if (clockwise_ordering)
47 std::reverse(indices.begin(), indices.end());
48 }
49 else if (reference_cell_type == ReferenceCells::Triangle)
50 {
51 indices = {deal2cgal.at(face->vertex_index(0)),
52 deal2cgal.at(face->vertex_index(1)),
53 deal2cgal.at(face->vertex_index(2))};
54
55 if (clockwise_ordering)
56 std::reverse(indices.begin(), indices.end());
57 }
58
59 else if (reference_cell_type == ReferenceCells::Quadrilateral)
60 {
61 if (clockwise_ordering)
62 {
63 indices = {deal2cgal.at(face->vertex_index(0)),
64 deal2cgal.at(face->vertex_index(2)),
65 deal2cgal.at(face->vertex_index(3)),
66 deal2cgal.at(face->vertex_index(1))};
67 }
68 else
69 {
70 indices = {deal2cgal.at(face->vertex_index(0)),
71 deal2cgal.at(face->vertex_index(1)),
72 deal2cgal.at(face->vertex_index(3)),
73 deal2cgal.at(face->vertex_index(2))};
74 }
75 }
76 else
78
79 [[maybe_unused]] const auto new_face = mesh.add_face(indices);
80 Assert(new_face != mesh.null_face(),
81 ExcInternalError("While trying to build a CGAL facet, "
82 "CGAL encountered a orientation problem that it "
83 "was not able to solve."));
84 }
85
86
87
88 template <typename dealiiFace, typename CGAL_Mesh>
89 void
90 map_vertices(
91 const dealiiFace &cell,
92 std::map<unsigned int, typename CGAL_Mesh::Vertex_index> &deal2cgal,
93 CGAL_Mesh &mesh)
94 {
95 for (const auto i : cell->vertex_indices())
96 {
97 if (deal2cgal.find(cell->vertex_index(i)) == deal2cgal.end())
98 {
99 deal2cgal[cell->vertex_index(i)] =
100 mesh.add_vertex(CGALWrappers::dealii_point_to_cgal_point<
101 typename CGAL_Mesh::Point>(cell->vertex(i)));
102 }
103 }
104 }
105} // namespace
106
107
108
109# ifndef DOXYGEN
110// Template implementations
111namespace CGALWrappers
112{
113 template <typename CGALPointType, int dim, int spacedim>
114 void
117 const Mapping<dim, spacedim> &mapping,
118 CGAL::Surface_mesh<CGALPointType> &mesh)
119 {
120 Assert(dim > 1, ExcImpossibleInDim(dim));
121 using Mesh = CGAL::Surface_mesh<CGALPointType>;
122 const auto &vertices = mapping.get_vertices(cell);
123 std::map<unsigned int, typename Mesh::Vertex_index> deal2cgal;
124
125 // Add all vertices to the mesh
126 // Store CGAL ordering
127 for (const auto &i : cell->vertex_indices())
128 deal2cgal[cell->vertex_index(i)] = mesh.add_vertex(
129 CGALWrappers::dealii_point_to_cgal_point<CGALPointType>(vertices[i]));
130
131 // Add faces
132 if (dim < 3)
133 // simplices and quads are allowable faces for CGAL
134 add_facet(cell, deal2cgal, mesh);
135 else
136 // in 3d, we build a surface mesh containing all the faces of the 3d cell.
137 // Simplices, Tetrahedrons, and Pyramids have their faces numbered in the
138 // same way as CGAL does (all faces of a bounding polyhedron are numbered
139 // counter-clockwise, so that their normal points outwards). Hexahedrons
140 // in deal.II have their faces numbered lexicographically, and one cannot
141 // deduce the direction of the normals by just looking at the vertices.
142 //
143 // In order for CGAL to be able to produce the right orientation, we need
144 // to reverse the order of the vertices for faces with even index.
145 // However, in order to allow for all kinds of meshes in 3d, including
146 // Moebius-loops, a deal.II face might even be rotated looking from one
147 // cell, whereas it is according to the standard when looking at it from
148 // the neighboring cell sharing that particular face. Therefore, when
149 // building a cgal face we must also take into account the fact that a
150 // face may have a non-standard orientation.
151 for (const auto &f : cell->face_indices())
152 {
153 // Check for standard orientation of faces
154 bool face_is_clockwise_oriented =
156 (f % 2 == 0);
157 // Make sure that we revert the orientation if required
158 if (cell->face_orientation(f) == false)
159 face_is_clockwise_oriented = !face_is_clockwise_oriented;
160 add_facet(cell->face(f), deal2cgal, mesh, face_is_clockwise_oriented);
161 }
162 }
163
164
165
166 template <typename CGALPointType, int dim, int spacedim>
167 void
169 const ::Triangulation<dim, spacedim> &tria,
170 CGAL::Surface_mesh<CGALPointType> &mesh)
171 {
172 Assert(tria.n_cells() > 0,
174 "Triangulation cannot be empty upon calling this function."));
175 Assert(mesh.is_empty(),
177 "The surface mesh must be empty upon calling this function."));
178
179 Assert(dim > 1, ExcImpossibleInDim(dim));
180 using Mesh = CGAL::Surface_mesh<CGALPointType>;
181 using Vertex_index = typename Mesh::Vertex_index;
182
183 std::map<unsigned int, Vertex_index> deal2cgal;
184 if constexpr (dim == 2)
185 {
186 for (const auto &cell : tria.active_cell_iterators())
187 {
188 map_vertices(cell, deal2cgal, mesh);
189 add_facet(cell, deal2cgal, mesh);
190 }
191 }
192 else if constexpr (dim == 3 && spacedim == 3)
193 {
194 for (const auto &cell : tria.active_cell_iterators())
195 {
196 for (const auto &f : cell->face_indices())
197 {
198 if (cell->face(f)->at_boundary())
199 {
200 map_vertices(cell->face(f), deal2cgal, mesh);
201 add_facet(cell->face(f),
202 deal2cgal,
203 mesh,
204 (f % 2 == 0 || cell->n_vertices() != 8));
205 }
206 }
207 }
208 }
209 else
210 {
211 Assert(false, ExcImpossibleInDimSpacedim(dim, spacedim));
212 }
213 }
214
215// Explicit instantiations.
216//
217// We don't build the instantiations.inst file if deal.II isn't
218// configured with CGAL, but doxygen doesn't know that and tries to
219// find that file anyway for parsing -- which then of course it fails
220// on. So exclude the following from doxygen consideration.
221# ifndef DOXYGEN
222# include "cgal/surface_mesh.inst"
223# endif
224
225} // namespace CGALWrappers
226# endif
227
228
230
231#endif
TriaIterator< TriaAccessor< dim - 1, dim, spacedim > > face(const unsigned int i) const
virtual boost::container::small_vector< Point< spacedim >, ReferenceCells::max_n_vertices< dim >() > get_vertices(const typename Triangulation< dim, spacedim >::cell_iterator &cell) const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > vertex_indices() const
unsigned int n_vertices() const
std_cxx20::ranges::iota_view< unsigned int, unsigned int > face_indices() const
unsigned int vertex_index(const unsigned int i) const
ReferenceCell reference_cell() const
#define DEAL_II_NAMESPACE_OPEN
Definition config.h:40
#define DEAL_II_NAMESPACE_CLOSE
Definition config.h:41
#define DEAL_II_ASSERT_UNREACHABLE()
static ::ExceptionBase & ExcImpossibleInDimSpacedim(int arg1, int arg2)
#define Assert(cond, exc)
static ::ExceptionBase & ExcImpossibleInDim(int arg1)
static ::ExceptionBase & ExcInternalError()
static ::ExceptionBase & ExcMessage(std::string arg1)
TriaIterator< CellAccessor< dim, spacedim > > cell_iterator
Definition tria.h:1557
bool face_orientation(const unsigned int face) const
MappingQ< dim, spacedim > StaticMappingQ1< dim, spacedim >::mapping
Definition mapping_q1.h:104
void dealii_tria_to_cgal_surface_mesh(const ::Triangulation< dim, spacedim > &triangulation, CGAL::Surface_mesh< CGALPointType > &mesh)
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)
constexpr ReferenceCell Triangle
constexpr ReferenceCell Hexahedron
constexpr ReferenceCell Quadrilateral
constexpr ReferenceCell Line