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

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

Definition

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

An angle-based regularization type for 2D segments that reinforces parallelism and orthogonality relationships.

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 input range 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_angles.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>
 Angle_regularization_2 (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. More...
 

Access

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

Groups

template<typename OutIterator >
OutIterator parallel_groups (OutIterator groups) const
 creates groups of indices, where each group represents parallel segments. More...
 
template<typename OutIterator >
OutIterator orthogonal_groups (OutIterator groups) const
 creates groups of indices, where each group represents orthogonal segments. More...
 

Miscellaneous

std::size_t number_of_modified_segments () const
 returns the number of modified segments.
 

Memory Management

void clear ()
 clears all internal data structures.
 

Constructor & Destructor Documentation

◆ Angle_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::Angle_regularization_2< GeomTraits, InputRange, SegmentMap >::Angle_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 angle deviation in degrees of a segment from its initial orientation
  • Type: GeomTraits::FT
  • Default: 25 degrees
  • an instance of SegmentMap that maps an item from input_range to GeomTraits::Segment_2
  • Default: SegmentMap()
Precondition
input_range.size() >= 2
maximum_angle >= 0 && maximum_angle <= 90

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::Angle_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

◆ orthogonal_groups()

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

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

This method calls Segments::orthogonal_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

◆ parallel_groups()

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

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

This method calls Segments::parallel_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::Angle_regularization_2< GeomTraits, InputRange, SegmentMap >::target ( const std::size_t  i,
const std::size_t  j 
) const

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

The target value is the angle difference between initial orientations of two segments i and j.

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

◆ update()

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

applies new orientations 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 angle that is added to the initial segment orientation.

Parameters
solutiona vector with angles in degrees
Precondition
solution.size() >= 1