CGAL 5.6.1 - Shape Regularization
CGAL::Shape_regularization::Segments::Delaunay_neighbor_query_2< GeomTraits, InputRange, SegmentMap > Class Template Reference

#include <CGAL/Shape_regularization/Segments/Delaunay_neighbor_query_2.h>

Definition

template<typename GeomTraits, typename InputRange, typename SegmentMap = CGAL::Identity_property_map<typename GeomTraits::Segment_2>>
class CGAL::Shape_regularization::Segments::Delaunay_neighbor_query_2< GeomTraits, InputRange, SegmentMap >

A neighbor query based on a Delaunay triangulation, which enables to find the nearest neighbors in a set of GeomTraits::Segment_2.

This class finds indices of the nearest neighbors of a query segment in a set of 2D segments.

The class first creates a Delaunay triangulation whose vertices are center points of the input segments and then enables to return for each query segment indices of the segments whose center points form the one-ring neighborhood of the corresponding query vertex in the triangulation.

Template Parameters
GeomTraitsa model of Kernel
InputRangea model of ConstRange whose iterator type is RandomAccessIterator
SegmentMapa model of ReadablePropertyMap whose key type is the value type of the InputRange and value type is GeomTraits::Segment_2. The default is CGAL::Identity_property_map<typename GeomTraits::Segment_2>.
Is Model Of:
NeighborQuery
Examples:
Shape_regularization/regularize_100_segments_angles.cpp, Shape_regularization/regularize_100_segments_offsets.cpp, Shape_regularization/regularize_15_segments.cpp, and Shape_regularization/regularize_real_data_2.cpp.

Initialization

template<typename NamedParameters = parameters::Default_named_parameters>
 Delaunay_neighbor_query_2 (const InputRange &input_range, const NamedParameters &np=parameters::default_values())
 initializes all internal data structures. More...
 
template<typename IndexRange >
void add_group (const IndexRange &index_range)
 inserts a group of segments from input_range and finds their neighbors within the group. More...
 

Access

void operator() (const std::size_t query_index, std::vector< std::size_t > &neighbors) const
 fills in a vector with indices of segments, which are direct neighbors of the query segment. More...
 

Memory Management

void clear ()
 clears all internal data structures.
 

Constructor & Destructor Documentation

◆ Delaunay_neighbor_query_2()

template<typename GeomTraits , typename InputRange , typename SegmentMap = CGAL::Identity_property_map<typename GeomTraits::Segment_2>>
template<typename NamedParameters = parameters::Default_named_parameters>
CGAL::Shape_regularization::Segments::Delaunay_neighbor_query_2< GeomTraits, InputRange, SegmentMap >::Delaunay_neighbor_query_2 ( const InputRange &  input_range,
const NamedParameters &  np = parameters::default_values() 
)

initializes all internal data structures.

Template Parameters
NamedParametersa sequence of Named Parameters
Parameters
input_rangea const range of 2D segments
npan optional sequence of Named Parameters among the ones listed below; this parameter can be omitted, the default values are then used
Optional Named Parameters
  • an instance of SegmentMap that maps an item from input_range to GeomTraits::Segment_2
  • Default: SegmentMap()
Precondition
input_range.size() >= 2

Member Function Documentation

◆ add_group()

template<typename GeomTraits , typename InputRange , typename SegmentMap = CGAL::Identity_property_map<typename GeomTraits::Segment_2>>
template<typename IndexRange >
void CGAL::Shape_regularization::Segments::Delaunay_neighbor_query_2< GeomTraits, InputRange, SegmentMap >::add_group ( const IndexRange &  index_range)

inserts a group of segments from input_range and finds their neighbors within the group.

Each group of segments is provided as a collection of their indices and only neighbors of segments within the group are computed that is no relationships between segments from different groups are taken into account.

The user must not use this method until he has meaningful groups of segments (see more in the user manual). By default, all segments are inserted as a group.

Template Parameters
IndexRangea model of ConstRange whose value type is std::size_t
Parameters
index_rangea const range of segment indices
Precondition
index_range.size() >= 2

◆ operator()()

template<typename GeomTraits , typename InputRange , typename SegmentMap = CGAL::Identity_property_map<typename GeomTraits::Segment_2>>
void CGAL::Shape_regularization::Segments::Delaunay_neighbor_query_2< GeomTraits, InputRange, SegmentMap >::operator() ( const std::size_t  query_index,
std::vector< std::size_t > &  neighbors 
) const

fills in a vector with indices of segments, which are direct neighbors of the query segment.

Parameters
query_indexan index of the query segment
neighborsindices of segments, which are direct neighbors of the query segment
Precondition
query_index >= 0 && query_index < input_range.size()