CGAL 5.3.1 - 2D Conforming Triangulations and Meshes

The main function to generate a mesh is refine_Delaunay_mesh_2().

The function lloyd_optimize_mesh_2() allows to optimize an existing mesh.

Functions

template<class CDT , class Criteria >
void CGAL::refine_Delaunay_mesh_2 (CDT &t, const Criteria &criteria=Criteria())
 Refines the default domain defined by a constrained Delaunay triangulation without seeds into a mesh satisfying the criteria defined by the traits criteria. More...
 
template<class CDT , class Criteria , class InputIterator >
void CGAL::refine_Delaunay_mesh_2 (CDT &t, InputIterator begin, InputIterator end, const Criteria &criteria=Criteria(), bool mark=false)
 Refines the default domain defined by a constrained Delaunay triangulation into a mesh satisfying the criteria defined by the traits criteria.The sequence [begin, end) gives a set of seeds points, that defines the domain to be meshed as follows. More...
 
template<typename CDT , typename PointIterator >
CGAL::Mesh_optimization_return_code CGAL::lloyd_optimize_mesh_2 (CDT &cdt, double parameters::time_limit=0, std::size_t parameters::max_iteration_number=0, double parameters::convergence=0.001, double parameters::freeze_bound=0.001, PointIterator parameters::seeds_begin=PointIterator(), PointIterator parameters::seeds_end=PointIterator(), bool parameters::mark=false)
 The function lloyd_optimize_mesh_2() is a mesh optimization process based on the minimization of a global energy function. More...
 
template<class CDT >
void CGAL::make_conforming_Delaunay_2 (CDT &t)
 Refines the constrained Delaunay triangulation t into a conforming Delaunay triangulation. More...
 
template<class CDT >
void CGAL::make_conforming_Gabriel_2 (CDT &t)
 Refines the constrained Delaunay triangulation t into a conforming Gabriel triangulation. More...
 

Function Documentation

◆ lloyd_optimize_mesh_2()

template<typename CDT , typename PointIterator >
CGAL::Mesh_optimization_return_code CGAL::lloyd_optimize_mesh_2 ( CDT &  cdt,
double parameters::time_limit  = 0,
std::size_t parameters::max_iteration_number  = 0,
double parameters::convergence  = 0.001,
double parameters::freeze_bound  = 0.001,
PointIterator parameters::seeds_begin  = PointIterator(),
PointIterator parameters::seeds_end  = PointIterator(),
bool parameters::mark  = false 
)

#include <CGAL/lloyd_optimize_mesh_2.h>

The function lloyd_optimize_mesh_2() is a mesh optimization process based on the minimization of a global energy function.

In lloyd_optimize_mesh_2(), the minimized global energy may be interpreted as the \( L^1\)-norm of the error achieved when the function \( x^2\) is interpolated on the mesh domain using a piecewise linear function which is linear in each cell of the Voronoi diagram of the mesh vertices.

The optimizer lloyd_optimize_mesh_2() works in iterative steps. At each iteration, mesh vertices are moved into positions that bring to zero the energy gradient and the Delaunay triangulation is updated. Vertices on the mesh boundaries are not moved.

Template Parameters
CDTis required to be or derive from CGAL::Constrained_Delaunay_triangulation_2, with vertex base and face base of its underlying TriangulationDataStructure_2 being models of DelaunayMeshVertexBase_2 and DelaunayMeshFaceBase_2, respectively. The argument cdt, passed by reference, provides the initial mesh and is modified by the algorithm to represent the final optimized mesh.
PointIteratormust be an iterator with value type Kernel::Point_2

The function has several optional parameters which are named parameters (we use the Boost.Parameter library). Therefore, when calling the function, the parameters can be provided in any order provided that the names of the parameters are used (see example at the bottom of this page).

Named Parameters

  • parameters::time_limit is used to set up, in seconds, a CPU time limit after which the optimization process is stopped. This time is measured using CGAL::Timer. The default value is 0 and means that there is no time limit.
    Precondition
    time_limit \( \geq\) 0
  • parameters::max_iteration_number sets a limit on the number of performed iterations. The default value of 0 means that there is no limit on the number of performed iterations.
    Precondition
    max_iteration_number \( \geq\) 0
  • parameters::convergence is a stopping criterion based on convergence: the optimization process is stopped, when at the last iteration, the displacement of any vertex is less than a given fraction of the length of the shortest edge incident to that vertex. The parameter convergence gives the threshold ratio.
    Precondition
    0 \( \leq\) convergence \( \leq\) 1
  • parameters::freeze_bound is designed to reduce running time of each optimization iteration. Any vertex that has a displacement less than a given fraction of the length of its shortest incident edge, is frozen (i.e. is not relocated). The parameter freeze_bound gives the threshold ratio. The default value is 0.001. If it is set to 0, freezing of vertices is disabled.
    Precondition
    0 \( \leq\) freeze_bound \( \leq\) 1
  • parameters::seeds_begin and parameters::seeds_end are begin and end input iterators to iterate on seed points. The sequence [parameters::seeds_begin, parameters::seeds_end) defines the domain in which the mesh was generated, and should be optimized.
  • parameters::mark. If mark is set to true, the mesh domain is the union of the bounded connected components including at least one seed. If mark is false, the domain is the union of the bounded components including no seed. Note that the unbounded component of the plane is never optimized. The default value is false.
Returns
The function lloyd_optimize_mesh_2() returns a value of type CGAL::Mesh_optimization_return_code which is:

Example

// Lloyd-smoothing until convergence reaches 0.01, freezing vertices which
// move less than 0.001*shortest_incident_edge_length
parameters::convergence=0.01,
parameters::freeze_bound=0.001);
See also
CGAL::Mesh_optimization_return_code
CGAL::refine_Delaunay_mesh_2()
Examples:
Mesh_2/mesh_optimization.cpp.

◆ make_conforming_Delaunay_2()

template<class CDT >
void CGAL::make_conforming_Delaunay_2 ( CDT &  t)

#include <CGAL/Triangulation_conformer_2.h>

Refines the constrained Delaunay triangulation t into a conforming Delaunay triangulation.

After a call to this function, all edges of t are Delaunay edges.

Template Parameters
CDTmust be a 2D constrained Delaunay triangulation and its geometric traits class must be a model of ConformingDelaunayTriangulationTraits_2.
Examples:
Mesh_2/conforming.cpp.

◆ make_conforming_Gabriel_2()

template<class CDT >
void CGAL::make_conforming_Gabriel_2 ( CDT &  t)

#include <CGAL/Triangulation_conformer_2.h>

Refines the constrained Delaunay triangulation t into a conforming Gabriel triangulation.

After a call to this function, all constrained edges of t have the Gabriel property: the circle that has \( e\) as diameter does not contain any vertex from the triangulation.

Template Parameters
CDTmust be a 2D constrained Delaunay triangulation and its geometric traits class must be a model of ConformingDelaunayTriangulationTraits_2.
Examples:
Mesh_2/conforming.cpp.

◆ refine_Delaunay_mesh_2() [1/2]

template<class CDT , class Criteria >
void CGAL::refine_Delaunay_mesh_2 ( CDT &  t,
const Criteria &  criteria = Criteria() 
)

#include <CGAL/Delaunay_mesher_2.h>

Refines the default domain defined by a constrained Delaunay triangulation without seeds into a mesh satisfying the criteria defined by the traits criteria.

The domain of the mesh covers all the connected components of the plane defined by the constrained edges of t, except for the unbounded component.

Template Parameters
CDTmust be 2D constrained Delaunay triangulation and its geometric traits class must be a model of DelaunayMeshTraits_2. The face of the constrained Delaunay triangulation must be a model of the concept DelaunayMeshFaceBase_2.
Criteriamust be a model of the concept MeshingCriteria_2. Criteria::Face_handle must be the same as CDT::Face_handle.
Examples:
Mesh_2/mesh_global.cpp, and Mesh_2/mesh_with_seeds.cpp.

◆ refine_Delaunay_mesh_2() [2/2]

template<class CDT , class Criteria , class InputIterator >
void CGAL::refine_Delaunay_mesh_2 ( CDT &  t,
InputIterator  begin,
InputIterator  end,
const Criteria &  criteria = Criteria(),
bool  mark = false 
)

#include <CGAL/Delaunay_mesher_2.h>

Refines the default domain defined by a constrained Delaunay triangulation into a mesh satisfying the criteria defined by the traits criteria.The sequence [begin, end) gives a set of seeds points, that defines the domain to be meshed as follows.

The constrained edges of t partition the plane into connected components. If mark==true, the mesh domain is the union of the bounded connected components including at least one seed. If mark==false, the domain is the union of the bounded components including no seed. Note that the unbounded component of the plane is never meshed.

Template Parameters
CDTmust be 2D constrained Delaunay triangulation and its geometric traits class must be a model of DelaunayMeshTraits_2. The face of the constrained Delaunay triangulation must be a model of the concept DelaunayMeshFaceBase_2.
Criteriamust be a model of the concept MeshingCriteria_2. Criteria::Face_handle must be the same as CDT::Face_handle.
InputIteratormust be an input iterator with value type CDT::Geom_traits::Point_2.