Function

CGAL::refine_Delaunay_mesh_2

template<class CDT, class Criteria>
void refine_Delaunay_mesh_2 ( CDT &t, 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. 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.
Precondition: The template parameter CDT must be a model of the concept ConstrainedDelaunayTriangulation_2. The geometric traits class of the constrained Delaunay triangulation must be a model of DelaunayMeshTraits_2.
Requirement: The face of the constrained Delaunay triangulation must be a model of the concept DelaunayMeshFaceBase_2. Criteria must be a model of the concept MeshingCriteria_2 and CDT::Face_handle must be the same as Criteria::Face_handle.

template <class CDT, class Criteria, class InputIterator>
void
refine_Delaunay_mesh_2 ( CDT& t,
InputIterator begin,
InputIterator end,
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. 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.
Requirement: The value_type of begin and end is CDT::Geom_traits::Point_2.