mirror of https://github.com/CGAL/cgal
338 lines
15 KiB
Plaintext
338 lines
15 KiB
Plaintext
namespace CGAL {
|
|
/*!
|
|
|
|
\mainpage User Manual
|
|
\anchor Chapter_2D_Conforming_Triangulations_and_Meshes
|
|
\anchor userchapter2DMeshes
|
|
\cgalAutoToc
|
|
\author Laurent Rineau
|
|
|
|
This package implements Shewchuk's algorithm \cgalCite{s-mgdsa-00} to construct
|
|
conforming triangulations and 2D meshes. Conforming triangulations will be
|
|
described in Section \ref secMesh_2_conforming_triangulation and
|
|
meshes in Section \ref secMesh_2_meshes.
|
|
|
|
\section secMesh_2_conforming_triangulation Conforming Triangulations
|
|
|
|
\subsection secMesh_2_conforming_definitions Definitions
|
|
|
|
A triangulation is a <I>Delaunay triangulation</I> if the circumscribing
|
|
circle of any facet of the triangulation contains no vertex in its
|
|
interior. A <I>constrained</I> Delaunay triangulation is a constrained
|
|
triangulation which is as much Delaunay as possible. The circumscribing
|
|
circle of any facet of a constrained Delaunay triangulation contains in its
|
|
interior no data point <I>visible</I> from the facet.
|
|
|
|
An edge is said to be a <I>Delaunay edge</I> if it is inscribed in an empty
|
|
circle (containing no data point in its interior). This edge is said to be a
|
|
<I>Gabriel edge</I> if its diametrical circle is empty.
|
|
|
|
A constrained Delaunay triangulation is said to be a <I>conforming
|
|
Delaunay triangulation</I> if every constrained edge is a Delaunay edge.
|
|
Because any edge in a constrained Delaunay triangulation is either a
|
|
Delaunay edge or a constrained edge, a conforming Delaunay triangulation is
|
|
in fact a Delaunay triangulation. The only difference is that some of the
|
|
edges are marked as constrained edges.
|
|
|
|
A constrained Delaunay triangulation is said to be a <I>conforming
|
|
Gabriel triangulation</I> if every constrained edge is a Gabriel edge. The
|
|
Gabriel property is stronger than the Delaunay property and each Gabriel
|
|
edge is a Delaunay edge. Conforming Gabriel triangulations are thus also
|
|
conforming Delaunay triangulations.
|
|
|
|
Any constrained Delaunay triangulation can be refined into a
|
|
conforming Delaunay triangulation or into a conforming Gabriel
|
|
triangulation by adding vertices, called <I>Steiner vertices</I>, on
|
|
constrained edges until they are decomposed into subconstraints small enough
|
|
to be Delaunay or Gabriel edges.
|
|
|
|
\subsection secMesh_2_building_conforming Building Conforming Triangulations
|
|
|
|
Constrained Delaunay triangulations can be refined into
|
|
conforming triangulations by the two following global functions:
|
|
|
|
\code{.cpp}
|
|
template<class CDT>
|
|
void make_conforming_Delaunay_2 (CDT& t)
|
|
\endcode
|
|
|
|
\code{.cpp}
|
|
template<class CDT>
|
|
void make_conforming_Gabriel_2 (CDT& t)
|
|
\endcode
|
|
|
|
|
|
In both cases, the template parameter `CDT` must be instantiated by a
|
|
constrained Delaunay triangulation class (see Chapter \ref
|
|
Chapter_2D_Triangulations "2D Triangulations").
|
|
|
|
The geometric traits of the constrained
|
|
Delaunay triangulation used to instantiate the parameter `CDT` has to
|
|
be a model of the concept `ConformingDelaunayTriangulationTraits_2`.
|
|
|
|
The constrained Delaunay triangulation `t` is passed by reference
|
|
and is refined into a conforming Delaunay triangulation or into a
|
|
conforming Gabriel triangulation by adding vertices. The user is advised to
|
|
make a copy of the input triangulation in the case where the original
|
|
triangulation has to be preserved for other computations
|
|
|
|
The algorithm used by `make_conforming_Delaunay_2()` and
|
|
`make_conforming_Gabriel_2()` builds internal data structures that would be
|
|
computed twice if the two functions are called consecutively on the same
|
|
triangulation. In order to avoid these data to be constructed twice, the
|
|
advanced user can use the class `Triangulation_conformer_2<CDT>` to
|
|
refine a constrained Delaunay triangulation into a conforming Delaunay
|
|
triangulation and then into a conforming Gabriel triangulation. For
|
|
additional control of the refinement algorithm, this class also provides
|
|
separate functions to insert one Steiner point at a time.
|
|
|
|
\subsection secMesh_2_example_making_conforming Example: Making a Triangulation Conforming Delaunay and Then Conforming Gabriel
|
|
|
|
This example inserts several segments into a constrained Delaunay
|
|
triangulation, makes it conforming Delaunay, and then conforming
|
|
Gabriel. At each step, the number of vertices of the triangulation is
|
|
printed.
|
|
|
|
\cgalExample{Mesh_2/conforming.cpp}
|
|
|
|
See \cgalFigureRef{Conformexampleconform}
|
|
|
|
\cgalFigureBegin{Conformexampleconform,example-conform-Delaunay-Gabriel.png}
|
|
Initial triangulation and the corresponding Delaunay and Gabriel triangulation.
|
|
\cgalFigureEnd
|
|
|
|
\cgalFigureBegin{ConformexampleconformDelaunay,example-conform-Delaunay.png}
|
|
The corresponding conforming Delaunay triangulation.
|
|
\cgalFigureEnd
|
|
|
|
\section secMesh_2_meshes Meshes
|
|
|
|
\subsection secMesh_2_meshes_definition Definitions
|
|
|
|
A mesh is a partition of a given region into simplices whose shapes
|
|
and sizes satisfy several criteria.
|
|
|
|
The domain is the region that the user wants to mesh. It has to be
|
|
a bounded region of the plane. The domain is defined by a <I>planar
|
|
straight line graph</I>, <span class="textsc">Pslg</span> for short, which is a set of segments
|
|
such that two segments in the set are either disjoint or share an
|
|
endpoint. The segments of the <span class="textsc">Pslg</span> are constraints that will be
|
|
represented by a union of edges in the mesh. The <span class="textsc">Pslg</span> can also
|
|
contain isolated points that will appear as vertices of the mesh.
|
|
|
|
The segments of the <span class="textsc">Pslg</span> are either segments of the
|
|
boundary or internals constraints. The segments of the <span class="textsc">Pslg</span> have to
|
|
cover the boundary of the domain.
|
|
|
|
The <span class="textsc">Pslg</span> divides the plane into several connected components. By
|
|
default, the domain is the union of the bounded connected components.
|
|
|
|
See \cgalFigureRef{Domain} for an example of a domain
|
|
defined without using seed points, and a possible mesh of it.
|
|
|
|
\cgalFigureBegin{Domain,domain-domain-mesh.png}
|
|
A domain defined without seed points and the generated mesh.
|
|
\cgalFigureEnd
|
|
|
|
|
|
|
|
The user can override this default by providing a set of seed points. Either
|
|
seed points mark components to be meshed or they mark components not to be
|
|
meshed (holes).
|
|
|
|
See
|
|
\cgalFigureRef{Domainseeds} for another domain defined with the same <span class="textsc">Pslg</span> and two seed points used to define holes.
|
|
In the corresponding mesh these
|
|
two holes are triangulated but not meshed.
|
|
|
|
|
|
\cgalFigureBegin{Domainseeds,domain-seeds-domain-seeds-mesh.png}
|
|
A domain with two seeds points defining holes and the generated mesh.
|
|
\cgalFigureEnd
|
|
|
|
\subsection secMesh_2_criteria Shape and Size Criteria
|
|
|
|
The shape criterion for triangles is a lower bound \f$ B\f$ on the ratio
|
|
between the circumradius and the shortest edge length. Such a bound
|
|
implies a lower bound of \f$ \arcsin{\frac{1}{2B}}\f$ on the minimum angle
|
|
of the triangle and an upper bound of \f$ \pi - 2* \arcsin{\frac{1}{2B}}\f$
|
|
on the maximum angle. Unfortunately, the termination of the algorithm
|
|
is guaranteed only if \f$ B \ge \sqrt{2}\f$, which corresponds to a lower
|
|
bound of \f$ 20.7\f$ degrees over the angles.
|
|
|
|
The size criterion can be any criterion that tends to prefer small
|
|
triangles. For example, the size criterion can be an upper bound on the
|
|
length of longest edge of triangles, or an upper bound on the radius of the
|
|
circumcircle. The size bound can vary over the domain. For example,
|
|
the size criterion could impose a small size for the triangles intersecting
|
|
a given line.
|
|
|
|
Both types of criteria are defined in an object `criteria` passed as
|
|
parameter of the meshing functions.
|
|
|
|
\subsection Mesh_2TheMeshingAlgorithm The Meshing Algorithm
|
|
|
|
The input to a meshing problem is a <span class="textsc">Pslg</span> and a set of seeds
|
|
describing the domain to be meshed, and a set of size and shape
|
|
criteria. The algorithm implemented in this package starts with a
|
|
constrained Delaunay triangulation of the input <span class="textsc">Pslg</span> and produces a
|
|
mesh using the Delaunay refinement method. This method inserts new vertices to
|
|
the triangulation, as far as possible from other vertices, and stops when the
|
|
criteria are satisfied.
|
|
|
|
If all angles between incident segments of the input <span class="textsc">Pslg</span>
|
|
are greater than \f$ 60\f$ degrees and if the bound on the
|
|
circumradius/edge ratio is greater than \f$ \sqrt{2}\f$,
|
|
the algorithm is guaranteed to terminate with a mesh
|
|
satisfying the size and shape criteria.
|
|
|
|
If some input angles are smaller than \f$ 60\f$ degrees, the algorithm will
|
|
end up with a mesh in which some triangles violate the criteria near small
|
|
input angles. This is unavoidable since small angles formed
|
|
by input segments cannot be suppressed. Furthermore, it has been
|
|
shown (\cgalCite{s-mgdsa-00}), that some domains with small input angles
|
|
cannot be meshed with angles even smaller than the small input angles.
|
|
Note that if the domain is a polygonal region, the resulting mesh will
|
|
satisfy size and shape criteria except for the small input angles.
|
|
In addition, the algorithm may succeed in producing meshes with a lower
|
|
angle bound greater than \f$ 20.7\f$ degrees, but there is no such guarantee.
|
|
|
|
\subsection secMesh_2_building_meshes Building Meshes
|
|
|
|
Meshes are obtained from
|
|
constrained Delaunay triangulations by calling the global function :
|
|
|
|
\code{.cpp}
|
|
template<class CDT, class Criteria>
|
|
void refine_Delaunay_mesh_2 (CDT &t, const Criteria& criteria)
|
|
\endcode
|
|
|
|
The template parameter `CDT` must be instantiated by a constrained
|
|
Delaunay triangulation class. In order to override the domain,
|
|
a version of this function has two more arguments that define a sequence of
|
|
seed points.
|
|
|
|
The geometric traits class of `CDT` has to be a
|
|
model of the concept `DelaunayMeshTraits_2`. This concept
|
|
refines the concept `ConformingDelaunayTriangulationTraits_2`
|
|
adding the geometric predicates and constructors. The template parameter
|
|
`Criteria` must be a model of `MeshingCriteria_2`. This concept
|
|
defines criteria that the triangles have to satisfy.
|
|
\cgal provides two models for this concept:
|
|
<UL>
|
|
<LI>`Delaunay_mesh_criteria_2<CDT>`, that defines a shape criterion
|
|
that bounds the minimum angle of triangles,
|
|
<LI>`Delaunay_mesh_size_criteria_2<CDT>`, that adds to the previous
|
|
criterion a bound on the maximum edge length.
|
|
</UL>
|
|
|
|
If the function `refine_Delaunay_mesh_2()` is called several times on the
|
|
same triangulation with different criteria, the algorithm rebuilds the
|
|
internal data structure used for meshing at every call. In order to avoid
|
|
rebuild the data structure at every call, the advanced user can
|
|
use the class `Delaunay_mesher_2<CDT>`. This class provides also step
|
|
by step functions. Those functions insert one vertex at a time.
|
|
|
|
Any object of type `Delaunay_mesher_2<CDT>` is constructed from a
|
|
reference to a `CDT`, and has several member functions to define the
|
|
domain to be meshed and to mesh the `CDT`. See the example given below
|
|
and the reference manual for details. Note that the `CDT` should not be
|
|
externally modified during the life time of the `Delaunay_mesher_2<CDT>`
|
|
object.
|
|
|
|
Once the mesh is constructed, one can determine which faces of the
|
|
triangulation are in the mesh domain using the `is_in_domain()` member
|
|
function of the face type (see the concept `DelaunayMeshFaceBase_2`).
|
|
|
|
|
|
\subsection secMesh_2_optimization Optimization of Meshes with Lloyd
|
|
|
|
The package also provides a global function that runs Lloyd optimization iterations on the
|
|
mesh generated by Delaunay refinement. The goal of this mesh optimization is to
|
|
improve the angles inside the mesh, and make them as close as possible to 60 degrees.
|
|
|
|
\code{.cpp}
|
|
|
|
template< class CDT >
|
|
Mesh_optimization_return_code lloyd_optimize_mesh_2(CDT& cdt);
|
|
|
|
\endcode
|
|
Note that this global function has more parameters (see details in reference pages) to tune the optimization process.
|
|
|
|
This optimization process alternates relocating vertices to the center of mass
|
|
of their Voronoi cells, and updating the Delaunay connectivity of the triangulation.
|
|
The center of mass is computed with respect to a sizing function that was designed to
|
|
preserve the local density of points in the mesh generated by Delaunay refinement.
|
|
|
|
See Figure \cgalFigureRef{Lloydfigure} for a mesh generated by `refine_Delaunay_mesh_2()` and optimized
|
|
with `lloyd_optimize_mesh_2()`. Figure \cgalFigureRef{LloydHistogramfigure}
|
|
shows the histogram of angles inside these meshes.
|
|
|
|
\cgalFigureBegin{Lloydfigure,lloyd.png}
|
|
(Left) Mesh generated by `refine_Delaunay_mesh_2()` for a uniform sizing criterion.
|
|
(Right) shows the same mesh after 100 iterations of Lloyd optimization.
|
|
\cgalFigureEnd
|
|
|
|
\cgalFigureBegin{LloydHistogramfigure,lloyd-histograms.png}
|
|
Histograms of angles inside the mesh after Delaunay refinement, and after 10 and 100 iterations
|
|
of Lloyd optimization.
|
|
After Delaunay refinement, angles are in the interval [28.5; 121.9] degrees.
|
|
After 10 iterations of Lloyd optimization, they are in [29.1; 110.8]. 100 iterations take them to [29.3; 109.9].
|
|
\cgalFigureEnd
|
|
|
|
|
|
\subsection Mesh_2Examples Examples
|
|
|
|
\subsubsection Mesh_2ExampleUsingtheGlobalFunction Example Using the Global Function
|
|
|
|
The following example inserts several segments into a constrained
|
|
triangulation and then meshes it using the global function
|
|
`refine_Delaunay_mesh_2()`. The size and shape criteria are the default ones
|
|
provided by the criteria class `Delaunay_mesh_criteria_2<K>`. No seeds are
|
|
given, meaning that the mesh domain covers the whole plane except the
|
|
unbounded component.
|
|
|
|
\cgalExample{Mesh_2/mesh_global.cpp}
|
|
|
|
\subsubsection Mesh_2ExampleUsingtheClass Example Using the Class Delaunay_mesher_2<CDT>
|
|
|
|
This example uses the class `Delaunay_mesher_2<CDT>` and calls
|
|
the `refine_mesh()` member function twice, changing the size and
|
|
shape criteria in between. In such a case, using twice the global
|
|
function `refine_Delaunay_mesh_2()` would be less efficient,
|
|
because some internal structures needed by the algorithm would be
|
|
built twice.
|
|
|
|
\cgalExample{Mesh_2/mesh_class.cpp}
|
|
|
|
\subsubsection Mesh_2ExampleUsingSeeds Example Using Seeds
|
|
|
|
This example uses the global function `refine_Delaunay_mesh_2()` but
|
|
defines a domain by using one seed. The size and shape criteria are the
|
|
default ones provided by the criteria class
|
|
`Delaunay_mesh_criteria_2<K>`.
|
|
|
|
Once the mesh is constructed, the `is_in_domain()` member function of
|
|
faces is used to count them.
|
|
|
|
\cgalExample{Mesh_2/mesh_with_seeds.cpp}
|
|
|
|
|
|
\subsubsection Mesh_2ExampleLloyd Example Using the Lloyd optimizer
|
|
|
|
This example uses the global function `lloyd_optimize_mesh_2()`.
|
|
The mesh is generated using the function `refine_Delaunay_mesh_2()` of `CGAL::Delaunay_mesher_2`,
|
|
and is then optimized using `lloyd_optimize_mesh_2()`. The optimization will stop
|
|
after 10 (set by `max_iteration_number`) iterations of alternating vertex relocations
|
|
and Delaunay connectivity updates. More termination conditions can be used and are detailed
|
|
in the Reference Manual.
|
|
|
|
|
|
\cgalExample{Mesh_2/mesh_optimization.cpp}
|
|
|
|
|
|
|
|
*/
|
|
} /* namespace CGAL */
|
|
|