CGAL 5.6 - 2D Generalized Barycentric Coordinates
CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap > Class Template Reference

#include <CGAL/Barycentric_coordinates_2/Delaunay_domain_2.h>

Definition

template<typename VertexRange, typename GeomTraits, typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
class CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >

2D Delaunay domain restricted to a simple polygon.

This class implements a discretized domain restricted to a simple polygon. The interior part of the input polygon is triangulated and refined with respect to the user-specified shape size parameter. The final triangulation is a constrained Delaunay triangulation, where the constraints are the polygon edges.

Internally, the package 2D Conforming Triangulations and Meshes is used. See it for more details.

Template Parameters
VertexRangea model of ConstRange whose iterator type is RandomAccessIterator
GeomTraitsa model of BarycentricTraits_2
PointMapa model of ReadablePropertyMap whose key type is VertexRange::value_type and value type is Point_2. The default is CGAL::Identity_property_map.
Is Model Of:
DiscretizedDomain_2
Examples:
Barycentric_coordinates_2/harmonic_coordinates.cpp, Barycentric_coordinates_2/shape_deformation.cpp, and Barycentric_coordinates_2/terrain_height_modeling.cpp.

Types

typedef GeomTraits::FT FT
 Number type.
 
typedef GeomTraits::Point_2 Point_2
 Point type.
 

Initialization

 Delaunay_domain_2 (const VertexRange &polygon, const GeomTraits traits=GeomTraits(), const PointMap point_map=PointMap())
 initializes all internal data structures. More...
 
template<typename PointRange >
void create (const FT max_edge_length, const PointRange &seeds)
 creates a refined Delaunay triangulation restricted to the input polygon. More...
 

Access

template<typename OutIterator >
OutIterator barycenters (OutIterator b_begin) const
 computes barycenters of all generated triangles. More...
 
std::size_t number_of_vertices () const
 returns the number of triangulation vertices. More...
 
const Point_2vertex (const std::size_t query_index) const
 returns a const reference to the triangulation vertex with the index query_index. More...
 
bool is_on_boundary (const std::size_t query_index) const
 controls if the triangulation vertex with the index query_index is on the polygon boundary. More...
 
void operator() (const std::size_t query_index, std::vector< std::size_t > &neighbors) const
 returns the one-ring neighborhood of the triangulation vertex with the index query_index. More...
 
void locate (const Point_2 &query, std::vector< std::size_t > &triangle) const
 locates a triangle that contains a given query point. More...
 

Memory Management

void clear ()
 clears all internal data structures.
 
void release_memory ()
 releases all memory that is used internally.
 

Constructor & Destructor Documentation

◆ Delaunay_domain_2()

template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >::Delaunay_domain_2 ( const VertexRange &  polygon,
const GeomTraits  traits = GeomTraits(),
const PointMap  point_map = PointMap() 
)

initializes all internal data structures.

Parameters
polygonan instance of VertexRange with the vertices of a simple polygon
traitsa traits class with geometric objects, predicates, and constructions; the default initialization is provided
point_mapan instance of PointMap that maps a vertex from polygon to Point_2; the default initialization is provided
Precondition
polygon.size() >= 3
polygon is simple

Member Function Documentation

◆ barycenters()

template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
template<typename OutIterator >
OutIterator CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >::barycenters ( OutIterator  b_begin) const

computes barycenters of all generated triangles.

Template Parameters
OutIteratora model of OutputIterator that accepts points of type Point_2
Parameters
b_beginthe beginning of the destination range with the computed barycenters
Returns
an output iterator to the element in the destination range, one past the last barycenter stored

◆ create()

template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
template<typename PointRange >
void CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >::create ( const FT  max_edge_length,
const PointRange &  seeds 
)

creates a refined Delaunay triangulation restricted to the input polygon.

After the construction is completed, the first n vertices are the polygon vertices, the next m vertices are the vertices generated along the polygon boundary, the last k vertices are the vertices generated in the interior part of the polygon.

Template Parameters
PointRangea model of ConstRange whose value type is GeomTraits::Point_2
Parameters
max_edge_lengthan upper bound on the length of the longest edge; see Delaunay_mesh_size_criteria_2 for more details
seedscontains seed points indicating, which parts of the input polygon should be partitioned and subdivided

◆ is_on_boundary()

template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
bool CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >::is_on_boundary ( const std::size_t  query_index) const

controls if the triangulation vertex with the index query_index is on the polygon boundary.

This function implements DiscretizedDomain_2::is_on_boundary().

Parameters
query_indexan index of the query vertex
Precondition
query_index >= 0 && query_index < number_of_vertices()

◆ locate()

template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
void CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >::locate ( const Point_2 query,
std::vector< std::size_t > &  triangle 
) const

locates a triangle that contains a given query point.

If triangle is empty, the query point does not belong to the domain.

This function implements DiscretizedDomain_2::locate().

Parameters
querya query point
trianglestores indices of the vertices, which form a triangle, that contains the query point

◆ number_of_vertices()

template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
std::size_t CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >::number_of_vertices ( ) const

returns the number of triangulation vertices.

This function implements DiscretizedDomain_2::number_of_vertices().

◆ operator()()

template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
void CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >::operator() ( const std::size_t  query_index,
std::vector< std::size_t > &  neighbors 
) const

returns the one-ring neighborhood of the triangulation vertex with the index query_index.

This function implements DiscretizedDomain_2::operator()().

Parameters
query_indexan index of the query vertex
neighborsstores indices of the vertices, which from the one-ring neighborhood of the query vertex
Precondition
query_index >= 0 && query_index < number_of_vertices()

◆ vertex()

template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
const Point_2& CGAL::Barycentric_coordinates::Delaunay_domain_2< VertexRange, GeomTraits, PointMap >::vertex ( const std::size_t  query_index) const

returns a const reference to the triangulation vertex with the index query_index.

This function implements DiscretizedDomain_2::vertex().

Parameters
query_indexan index of the requested vertex
Precondition
query_index >= 0 && query_index < number_of_vertices()