wip round 6 Mael

This commit is contained in:
Jane Tournois 2025-06-13 16:17:36 +02:00
parent 2675d5c0c1
commit 9fb6dafd39
2 changed files with 35 additions and 31 deletions

View File

@ -15,10 +15,10 @@ namespace CGAL {
\section CT_3_CCDT_3 Constrained Triangulations in 3D
3D triangulations partition space and are useful in many applications. In some cases, it is
important to ensure that specific faces, such as sharp features, are preserved in the output.
important to ensure that specific faces, such as those representing sharp features, appear in the output.
When a triangulation exactly respects these constraints, it is called a _constrained_ triangulation.
However, it is sometimes only possible to preserve the geometry of the constraints, not their exact
combinatorics; in such cases, additional points—called Steiner points—must be inserted. This process
However, it is sometimes only possible to preserve the geometry of the constraints, but not their exact
combinatorics. In such cases, additional points, called _Steiner_ _points_, must be inserted. This process
results in a _conforming_ triangulation.
This package implements an algorithm for constructing conforming triangulations of 3D polygonal
@ -28,7 +28,8 @@ as described in the chapter \ref PkgTriangulation3.
The article by Cohen-Steiner et al. \cgalCite{cgal:cohen2002conforming} discusses the problem of
constructing conforming Delaunay triangulations and proposes an algorithm to address it.
Si's work \cgalCite{cgal:si2008cdt3} presents an algorithm for computing conforming constrained
Si's work \cgalCite{si2005meshing}, \cgalCite{cgal:si2008cdt3}, \cgalCite{si2015tetgen},
presents an algorithm for computing conforming constrained
Delaunay triangulations in 3D.
\section CT_3_definitions Definitions
@ -55,7 +56,7 @@ Polygons in a PLC may be non-convex, may have holes, and may have arbitrarily ma
<img src="plc.png" style="max-width:60%;"/>
</center>
\cgalFigureCaptionBegin{CT_3_plc_fig}
A Piecewise Linear Complex composed of planar faces connected by edges and vertices.
A Piecewise Linear Complex, composed of planar faces connected by edges and vertices.
\cgalFigureCaptionEnd
@ -94,8 +95,9 @@ version of the original PLC that admits a constrained Delaunay triangulation. Th
achieved by adding Steiner vertices to the input edges and polygons. The constrained triangulation
built on this refined PLC is known as a _conforming constrained Delaunay triangulation_ (CCDT for short).
The algorithm implemented in this package is based on the work of Hang Si \cgalCite{si2005meshing},
\cgalCite{si2015tetgen}. Steiner vertices are inserted on the input edges and polygons of the PLC to
The algorithm implemented in this package is based on the work of Hang Si
\cgalCite{si2005meshing}, \cgalCite{cgal:si2008cdt3}, \cgalCite{si2015tetgen}.
Steiner vertices are inserted on the input edges and polygons of the PLC to
ensure it can be tetrahedralized.
\cgalFigureRef{CT_3_plc2cdt_fig} illustrates an example of a conforming constrained
@ -138,13 +140,13 @@ This package also provides a way to group polygons into distinct surface patches
Each polygon can be assigned a _patch_ identifier, allowing multiple polygons to form a continuous surface patch,
which may include holes. When these patches are planar and satisfy the necessary geometric conditions,
they can be used to construct a conforming constrained Delaunay triangulation.
When a property map is provided:
- The vertices of the PLC are those from the original surface mesh or polygon soup.
When a property map is provided, the input PLC is then represented by its faces, edges and vertices as follows:
- Each face of the PLC is defined as the union of input polygons sharing the same patch identifier;
- The edges of the PLC are those belonging to the surface mesh or polygon soup that have only one
adjacent face, specifically those marking the boundaries of patches.
- The surface patches themselves serve as the polygons (faces) in the resulting representation.
adjacent face, specifically those marking the boundaries of patches;
- The vertices of the PLC are the ones lying at the boundaries of the surface patches in the original surface mesh or polygon soup.
\subsection CT_3_api_classes Classes
\subsection CT_3_api API
This package provides a primary class, `CGAL::Conforming_constrained_Delaunay_triangulation_3`.
This class is templated by a
@ -155,8 +157,6 @@ In addition to the main class, the package includes several auxiliary classes th
of vertices, cells, and associated metadata used within the triangulation. These supporting classes
enable users to extend or adapt the triangulation data structure to their specific needs.
\subsection CT_3_api_functions Functions
Two overloads of the constructor function \link
PkgConstrainedTriangulation3FunctionsPolygonSoupOrMesh
`CGAL::make_conforming_constrained_Delaunay_triangulation_3()`\endlink are provided to facilitate the creation
@ -201,6 +201,16 @@ When the named parameter `plc_face_id` is specified, each constrained facet in t
is assigned to the corresponding input PLC face, as identified by the provided property map.
If this parameter is not specified, each input polygon (or PLC face) is assigned a unique face index.
Figure \cgalFigureRef{CT_3_ccdt_examples_fig} shows the benefit of using the `plc_face_id` property map.
On the last line of the figure, the input PLC is enriched with a segmentation of the planar faces,
provided via the `plc_face_id` property map. In the resulting conforming constrained Delaunay triangulation,
only the boundary edges of the PLC faces are constrained, while the other edges never get inserted as
edges of the 3D triangulation.
Without the `plc_face_id` property map, all edges of the PLC faces are constrained,
each PLC face is considered as a constraint,
possibly resulting in a 3D triangulation with surfaces that are more refined than necessary.
\cgalExample{Constrained_triangulation_3/conforming_constrained_Delaunay_triangulation_3_fimap.cpp}
\cgalFigureRef{CT_3_ccdt_examples_fig} shows the input and output of this triangulation construction example.
@ -228,33 +238,26 @@ and does not contain self-intersections.
Several preprocessing functions are available in the \ref PkgPolygonMeshProcessing package to help ensure these
preconditions are met.
\subsubsection CT_3_example_ccdt_autorefinement Autorefinement of the Input Mesh
The following example shows how autorefinement can be used to preprocess the input mesh,
when it is self-intersecting.
\cgalExample{Constrained_triangulation_3/ccdt_3_after_autorefinement.cpp }
\subsubsection CT_3_example_self_intersecting_non_triangles Non-Triangulated Self-Intersecting Input Mesh
The following example demonstrates how to construct a conforming constrained Delaunay triangulation from
an input mesh that is not triangulated and may contain self-intersections.
an input mesh that is not triangulated and may contain self-intersections,
using autorefinement.
\cgalExample{Constrained_triangulation_3/ccdt_3_preprocessing.cpp }
The function
\link CGAL::Polygon_mesh_processing::does_self_intersect(const FaceRange&, const TriangleMesh&, const NamedParameters&) `CGAL::Polygon_mesh_processing::does_self_intersect()` \endlink
can be used to detect self-intersections,
is used to detect self-intersections,
but it requires the input mesh to be triangulated. Therefore, the input mesh must first be triangulated
using
\link CGAL::Polygon_mesh_processing::triangulate_faces(FaceRange,PolygonMesh&,const NamedParameters&) `CGAL::Polygon_mesh_processing::triangulate_faces()`\endlink
before performing the self-intersection check.
If self-intersections are found, the triangulated mesh is converted into a triangle soup, which can then be
If self-intersections are found, the triangulated mesh is converted into a triangle soup, which is then
processed with
\link CGAL::Polygon_mesh_processing::autorefine_triangle_soup(PointRange&,TriangleRange&,const NamedParameters&) `CGAL::Polygon_mesh_processing::autorefine_triangle_soup()`\endlink
to resolve the self-intersections.
\cgalExample{Constrained_triangulation_3/ccdt_3_preprocessing.cpp }
In case
\subsection CT_3_example_remesh Remeshing a Conforming Constrained Delaunay Triangulation
@ -283,8 +286,6 @@ The input data are:
<li> the input PLC, remeshed and not triangulated, with a segmentation of the surface patches, done with `CGAL::Polygon_mesh_processing::remesh_planar_patches()` as in example \ref CT_3_example_ccdt_fimap;</li>
<li> the input PLC, isotropically remeshed, with a segmentation of the surface patches, done with `CGAL::Polygon_mesh_processing::region_growing_of_planes_on_faces()` as in example \ref CT_3_example_ccdt_region_growing_fimap.</li>
</ul>
It is interesting to note that the conforming constrained Delaunay triangulations built from the same input PLC, but with different
preprocessing steps, can yield different results.
\cgalFigureAnchor{CT_3_ccdt_examples_fig}

View File

@ -43,8 +43,11 @@ int main(int argc, char* argv[])
std::cout << "Number of input triangles after autorefine: "
<< polygons.size() << "\n";
if(PMP::does_polygon_soup_self_intersect(points, polygons)) {
if(PMP::does_polygon_soup_self_intersect(points, polygons))
{
std::cerr << "Error: mesh still self-intersects after autorefine\n";
std::cerr << "Run autorefinement again with an exact kernel ";
std::cerr << "such as CGAL::Exact_predicates_exact_constructions_kernel \n";
return EXIT_FAILURE;
}