CGAL 5.6.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<typename CDT , typename NamedParameters = CGAL::parameters::Default_named_parameters>
void CGAL::refine_Delaunay_mesh_2 (CDT &t, const NamedParameters &np)
 refines the domain defined by a constrained Delaunay triangulation into a mesh satisfying the criteria defined by the traits criteria. More...
 
template<class CDT , class Criteria >
void CGAL::refine_Delaunay_mesh_2 (CDT &t, const Criteria &criteria=Criteria())
 
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)
 
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...
 
template<typename CDT , typename NamedParameters = CGAL::parameters::Default_named_parameters>
Mesh_optimization_return_code CGAL::lloyd_optimize_mesh_2 (CDT &cdt, const NamedParameters &np=parameters::default_values())
 The function lloyd_optimize_mesh_2() is a mesh optimization process based on the minimization of a global energy function. More...
 

Function Documentation

◆ lloyd_optimize_mesh_2()

template<typename CDT , typename NamedParameters = CGAL::parameters::Default_named_parameters>
Mesh_optimization_return_code CGAL::lloyd_optimize_mesh_2 ( CDT &  cdt,
const NamedParameters &  np = parameters::default_values() 
)

#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.
NamedParametersa sequence of Named Parameters
Parameters
cdtthe initial mesh that will be modified by the algorithm to represent the final optimized mesh.
npan optional sequence of Named Parameters among the ones listed below:
Optional Named Parameters
  • 2D points used to define the domain to mesh. If seeds_are_in_domain==true, the mesh domain is the union of the bounded connected components including at least one seed. If seeds_are_in_domain==false, the domain is the union of the bounded components including no seed.
  • Type: a class model of ConstRange whose iterator is a model of InputIterator with CDT::Point_2 as value type.
  • Default: No seed.

  • specified if seeds indicate bounded connected components inside or outside of the domain.
  • Type: bool
  • Default: false

  • limit on the number of performed iterations. 0 means that there is no limit on the number of performed iterations.
  • Precondition: max_iteration_number >=0
  • Type: int
  • Default: 0

  • CPU time limit (in seconds) after which the optimization process is stopped. This time is measured using CGAL::Real_timer. 0 means that there is no time limit.
  • Type: double
  • Precondition: time_limit \( \geq\) 0
  • Default: 0

  • 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. If it is set to 0, freezing of vertices is disabled.
  • Precondition: 0<= freeze_bound <=1
  • Type: double
  • Default: 0.001

  • threshold ratio of 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.
  • Precondition: 0 <=convergence <= 1
  • Type: double
  • Default: 0.001

Returns
an enum value providing some information about the outcome of the algorithm.
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/3]

template<typename CDT , typename NamedParameters = CGAL::parameters::Default_named_parameters>
void CGAL::refine_Delaunay_mesh_2 ( CDT &  t,
const NamedParameters &  np 
)

#include <CGAL/Delaunay_mesher_2.h>

refines the domain defined by a constrained Delaunay triangulation into a mesh satisfying the criteria defined by the traits criteria.

Note that the unbounded component of the plane is never meshed.

Template Parameters
CDTmust be a 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.
NamedParametersa sequence of Named Parameters
Parameters
tthe constrained Delaunay triangulation that will be refined
npan optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • 2D points used to define the domain to mesh. If seeds_are_in_domain==true, the mesh domain is the union of the bounded connected components including at least one seed. If seeds_are_in_domain==false, the domain is the union of the bounded components including no seed.
  • Type: a class model of ConstRange with iterator being a model of InputIterator with CDT::Point as value type.
  • Default: No seed.
  • If true the domain is considered as already initialized (see DelaunayMeshFaceBase_2::is_in_domain()). If false and no seeds are provided, 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.
  • Type: bool
  • Default: false
  • specified if seeds indicates bounded connected components inside or outside of the domain.
  • Type: bool
  • Default: false

Examples:
Mesh_2/mesh_global.cpp, Mesh_2/mesh_marked_domain.cpp, and Mesh_2/mesh_with_seeds.cpp.

◆ refine_Delaunay_mesh_2() [2/3]

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

#include <CGAL/Delaunay_mesher_2.h>

Deprecated:
This function is deprecated since CGAL 5.6

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 a 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.

◆ refine_Delaunay_mesh_2() [3/3]

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>

Deprecated:
This function is deprecated since CGAL 5.6

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.