CGAL 5.6 - Shape Regularization
Loading...
Searching...
No Matches
CGAL::Shape_regularization::Segments::Offset_regularization_2< GeomTraits, InputRange, SegmentMap > Class Template Reference

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

Definition

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

An offset-based regularization type for 2D segments that reinforces collinearity relationships.

The input groups of segments should each contain (near) parallel segments. In order to achieve that, one may use the class Segments::Angle_regularization_2 or the function Segments::parallel_groups(). Each group of parallel segments may be inserted using the method add_group().

Template Parameters
GeomTraitsa model of Kernel
InputRangea model of Range whose iterator type is RandomAccessIterator
SegmentMapa model of ReadWritePropertyMap 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:
RegularizationType
Examples
Shape_regularization/regularize_100_segments_offsets.cpp, Shape_regularization/regularize_15_segments.cpp, and Shape_regularization/regularize_real_data_2.cpp.

Types

typedef GeomTraits::FT FT
 Number type.
 

Initialization

template<typename NamedParameters = parameters::Default_named_parameters>
 Offset_regularization_2 (InputRange &input_range, const NamedParameters &np=parameters::default_values())
 initializes all internal data structures.
 
template<typename IndexRange >
void add_group (const IndexRange &index_range)
 inserts a group of segments from input_range.
 

Access

FT target (const std::size_t i, const std::size_t j) const
 calculates the target value between 2 parallel segments, which are direct neighbors to each other.
 
FT bound (const std::size_t) const
 returns maximum_offset.
 
void update (const std::vector< FT > &solution)
 applies new positions computed by the QuadraticProgramTraits to the initial segments.
 

Groups

template<typename OutIterator >
OutIterator collinear_groups (OutIterator groups) const
 creates groups of indices, where each group represents collinear segments.
 

Miscellaneous

std::size_t number_of_modified_segments () const
 returns the number of modified segments.
 
template<typename OutIterator >
OutIterator unique_segments (OutIterator segments) const
 returns segments, which best-fit collinear groups.
 

Memory Management

void clear ()
 clears all internal data structures.
 

Constructor & Destructor Documentation

◆ Offset_regularization_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::Offset_regularization_2< GeomTraits, InputRange, SegmentMap >::Offset_regularization_2 ( InputRange &  input_range,
const NamedParameters &  np = parameters::default_values() 
)

initializes all internal data structures.

Template Parameters
NamedParametersa sequence of Named Parameters
Parameters
input_rangea range of 2D segments to be regularized
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
  • maximum allowed orthogonal distance between two parallel segments such that they are considered to be collinear
  • Type: GeomTraits::FT
  • Default: 0.5 unit length
  • an instance of SegmentMap that maps an item from input_range to GeomTraits::Segment_2
  • Default: SegmentMap()
Precondition
input_range.size() >= 2
maximum_offset >= 0

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::Offset_regularization_2< GeomTraits, InputRange, SegmentMap >::add_group ( const IndexRange &  index_range)

inserts a group of segments from input_range.

Each group of segments is provided as a collection of their indices and only segments within the group are being regularized 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

◆ collinear_groups()

template<typename GeomTraits , typename InputRange , typename SegmentMap = CGAL::Identity_property_map<typename GeomTraits::Segment_2>>
template<typename OutIterator >
OutIterator CGAL::Shape_regularization::Segments::Offset_regularization_2< GeomTraits, InputRange, SegmentMap >::collinear_groups ( OutIterator  groups) const

creates groups of indices, where each group represents collinear segments.

This method calls Segments::collinear_groups().

Template Parameters
OutIteratora model of OutputIterator that accepts elements of type std::vector<std::size_t>
Parameters
groupsan output iterator with groups of segment indices

◆ target()

template<typename GeomTraits , typename InputRange , typename SegmentMap = CGAL::Identity_property_map<typename GeomTraits::Segment_2>>
FT CGAL::Shape_regularization::Segments::Offset_regularization_2< GeomTraits, InputRange, SegmentMap >::target ( const std::size_t  i,
const std::size_t  j 
) const

calculates the target value between 2 parallel segments, which are direct neighbors to each other.

The target value is the distance between two parallel segments i and j, where the distance is defined as the distance between the midpoint of the ith segment and the projection of this point onto the supporting line of the jth segment.

Parameters
iindex of the first segment
jindex of the second segment
Precondition
i >= 0 && i < input_range.size()
j >= 0 && j < input_range.size()

◆ unique_segments()

template<typename GeomTraits , typename InputRange , typename SegmentMap = CGAL::Identity_property_map<typename GeomTraits::Segment_2>>
template<typename OutIterator >
OutIterator CGAL::Shape_regularization::Segments::Offset_regularization_2< GeomTraits, InputRange, SegmentMap >::unique_segments ( OutIterator  segments) const

returns segments, which best-fit collinear groups.

This method first calls Segments::collinear_groups() and then substitutes each group of collinear segments by an average segment. The number of returned segments is the number of detected collinear groups.

Template Parameters
OutIteratora model of OutputIterator that accepts segments of type GeomTraits::Segment_2
Parameters
segmentsan output iterator with the simplified segments

◆ update()

template<typename GeomTraits , typename InputRange , typename SegmentMap = CGAL::Identity_property_map<typename GeomTraits::Segment_2>>
void CGAL::Shape_regularization::Segments::Offset_regularization_2< GeomTraits, InputRange, SegmentMap >::update ( const std::vector< FT > &  solution)

applies new positions computed by the QuadraticProgramTraits to the initial segments.

Number of values in solution equals to the number n of segments being regularized + the number m of neighbor pairs between these segments. Each of n values is an offset that is added to the initial segment position.

Parameters
solutiona vector with offsets in unit lengths
Precondition
solution.size() >= 1