#include <CGAL/Barycentric_coordinates_2/Delaunay_domain_2.h>
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
-
- 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.
|
typedef GeomTraits::FT | FT |
| Number type.
|
|
typedef GeomTraits::Point_2 | Point_2 |
| Point type.
|
|
|
| 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...
|
|
|
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_2 & | vertex (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...
|
|
|
void | clear () |
| clears all internal data structures.
|
|
void | release_memory () |
| releases all memory that is used internally.
|
|
◆ Delaunay_domain_2()
template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
initializes all internal data structures.
- Parameters
-
polygon | an instance of VertexRange with the vertices of a simple polygon |
traits | a traits class with geometric objects, predicates, and constructions; the default initialization is provided |
point_map | an instance of PointMap that maps a vertex from polygon to Point_2 ; the default initialization is provided |
- Precondition
- polygon.size() >= 3
-
polygon is simple
◆ barycenters()
template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
template<typename OutIterator >
computes barycenters of all generated triangles.
- Template Parameters
-
- Parameters
-
b_begin | the 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 >
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
-
PointRange | a model of ConstRange whose value type is GeomTraits::Point_2 |
- Parameters
-
max_edge_length | an upper bound on the length of the longest edge; see Delaunay_mesh_size_criteria_2 for more details |
seeds | contains 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>>
◆ locate()
template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
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
-
query | a query point |
triangle | stores 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>>
◆ operator()()
template<typename VertexRange , typename GeomTraits , typename PointMap = CGAL::Identity_property_map<typename GeomTraits::Point_2>>
returns the one-ring neighborhood of the triangulation vertex with the index query_index
.
This function implements DiscretizedDomain_2::operator()()
.
- Parameters
-
query_index | an index of the query vertex |
neighbors | stores 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>>
returns a const reference to the triangulation vertex with the index query_index
.
This function implements DiscretizedDomain_2::vertex()
.
- Parameters
-
query_index | an index of the requested vertex |
- Precondition
- query_index >= 0 && query_index < number_of_vertices()