CGAL 5.5.2 - Shape Regularization

Models and functions that can be used when regularizing segments.

Classes

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. More...
 
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. More...
 
class  CGAL::Shape_regularization::Segments::Offset_regularization_2< GeomTraits, InputRange, SegmentMap >
 An offset-based regularization type for 2D segments that reinforces collinearity relationships. More...
 

Functions

template<typename InputRange , typename NeighQuery , typename RegType , typename QPSolver , typename NamedParameters = parameters::Default_named_parameters>
void CGAL::Shape_regularization::Segments::regularize_segments (InputRange &input_range, NeighQuery &neighbor_query, RegType &regularization_type, QPSolver &quadratic_program, const NamedParameters &np=parameters::default_values())
 regularizes a set of 2D segments. More...
 
template<typename InputRange >
void CGAL::Shape_regularization::Segments::regularize_angles (InputRange &input_range)
 regularizes angles in a set of 2D segments. More...
 
template<typename InputRange >
void CGAL::Shape_regularization::Segments::regularize_offsets (InputRange &input_range)
 regularizes offsets in a set of 2D segments. More...
 
template<typename InputRange , typename OutIterator , typename NamedParameters = parameters::Default_named_parameters>
OutIterator CGAL::Shape_regularization::Segments::parallel_groups (const InputRange &input_range, OutIterator groups, const NamedParameters &np=parameters::default_values())
 finds groups of parallel segments in a set of 2D segments. More...
 
template<typename InputRange , typename OutIterator , typename NamedParameters = parameters::Default_named_parameters>
OutIterator CGAL::Shape_regularization::Segments::collinear_groups (const InputRange &input_range, OutIterator groups, const NamedParameters &np=parameters::default_values())
 finds groups of collinear segments in a set of 2D segments. More...
 
template<typename InputRange , typename OutIterator , typename NamedParameters = parameters::Default_named_parameters>
OutIterator CGAL::Shape_regularization::Segments::orthogonal_groups (const InputRange &input_range, OutIterator groups, const NamedParameters &np=parameters::default_values())
 finds groups of orthogonal segments in a set of 2D segments. More...
 
template<typename InputRange , typename OutIterator , typename NamedParameters = parameters::Default_named_parameters>
OutIterator CGAL::Shape_regularization::Segments::unique_segments (const InputRange &input_range, OutIterator segments, const NamedParameters &np=parameters::default_values())
 substitutes groups of 2D collinear segments by average segments. More...
 

Function Documentation

◆ collinear_groups()

template<typename InputRange , typename OutIterator , typename NamedParameters = parameters::Default_named_parameters>
OutIterator CGAL::Shape_regularization::Segments::collinear_groups ( const InputRange &  input_range,
OutIterator  groups,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Shape_regularization/regularize_segments.h>

finds groups of collinear segments in a set of 2D segments.

This function enables to find groups of near collinear segments in a set of 2D segments. The groups are returned as vectors of indices. This algorithm first finds the groups of parallel segments using the function Segments::parallel_groups() and then splits these groups into groups of collinear segments.

This function does not regularize input segments, but only groups them.

Template Parameters
InputRangea model of ConstRange whose iterator type is RandomAccessIterator
OutIteratora model of OutputIterator that accepts elements of type std::vector<std::size_t>
NamedParametersa sequence of Named Parameters
Parameters
input_rangea const range of input segments
groupsan output iterator with groups of segment indices
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.2 unit length
  • indicates whether the order of input segments should be preserved or not
  • Type: boolean
  • Default: false
  • a property map that maps an item from input_range to GeomTraits::Segment_2
  • Type: a model of ReadablePropertyMap whose key type is the value type of the input range and value type is GeomTraits::Segment_2
  • Default: CGAL::Identity_property_map
Returns
an output iterator to the element in the destination range, one past the last group stored
Precondition
input_range.size() >= 1
maximum_offset >= 0

◆ orthogonal_groups()

template<typename InputRange , typename OutIterator , typename NamedParameters = parameters::Default_named_parameters>
OutIterator CGAL::Shape_regularization::Segments::orthogonal_groups ( const InputRange &  input_range,
OutIterator  groups,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Shape_regularization/regularize_segments.h>

finds groups of orthogonal segments in a set of 2D segments.

This function enables to find groups of near orthogonal segments in a set of 2D segments. The groups are returned as vectors of indices. This algorithm first finds the groups of parallel segments using the function Segments::parallel_groups() and then merges these groups into groups of orthogonal segments.

This function does not regularize input segments, but only groups them.

Template Parameters
InputRangea model of ConstRange whose iterator type is RandomAccessIterator
OutIteratora model of OutputIterator that accepts elements of type std::vector<std::size_t>
NamedParametersa sequence of Named Parameters
Parameters
input_rangea const range of input segments
groupsan output iterator with groups of segment indices
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 between two segments such that they are considered to be parallel or orthogonal
  • Type: GeomTraits::FT
  • Default: 5 degrees
  • indicates whether the order of input segments should be preserved or not
  • Type: boolean
  • Default: false
  • a property map that maps an item from input_range to GeomTraits::Segment_2
  • Type: a model of ReadablePropertyMap whose key type is the value type of the input range and value type is GeomTraits::Segment_2
  • Default: CGAL::Identity_property_map
Returns
an output iterator to the element in the destination range, one past the last group stored
Precondition
input_range.size() >= 1
maximum_angle >= 0 && maximum_angle <= 90

◆ parallel_groups()

template<typename InputRange , typename OutIterator , typename NamedParameters = parameters::Default_named_parameters>
OutIterator CGAL::Shape_regularization::Segments::parallel_groups ( const InputRange &  input_range,
OutIterator  groups,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Shape_regularization/regularize_segments.h>

finds groups of parallel segments in a set of 2D segments.

This function enables to find groups of near parallel segments in a set of 2D segments. The groups are returned as vectors of indices. Note that two segments may be included at the same group even if they are far away from each other. This algorithm concerns only the angle relationship among segments, but not the distance.

This function does not regularize input segments, but only groups them.

Template Parameters
InputRangea model of ConstRange whose iterator type is RandomAccessIterator
OutIteratora model of OutputIterator that accepts elements of type std::vector<std::size_t>
NamedParametersa sequence of Named Parameters
Parameters
input_rangea const range of input segments
groupsan output iterator with groups of segment indices
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 between two segments such that they are considered to be parallel
  • Type: GeomTraits::FT
  • Default: 5 degrees
  • indicates whether the order of input segments should be preserved or not
  • Type: boolean
  • Default: false
  • a property map that maps an item from input_range to GeomTraits::Segment_2
  • Type: a model of ReadablePropertyMap whose key type is the value type of the input range and value type is GeomTraits::Segment_2
  • Default: CGAL::Identity_property_map
Returns
an output iterator to the element in the destination range, one past the last group stored
Precondition
input_range.size() >= 1
maximum_angle >= 0 && maximum_angle <= 90
Examples:
Shape_regularization/regularize_100_segments_offsets.cpp.

◆ regularize_angles()

template<typename InputRange >
void CGAL::Shape_regularization::Segments::regularize_angles ( InputRange &  input_range)

#include <CGAL/Shape_regularization/regularize_segments.h>

regularizes angles in a set of 2D segments.

Given a set of unordered 2D segments, this function enables to reinforce two types of regularities among these segments:

  • Parallelism: segments, which are detected as near parallel, are made exactly parallel.
  • Orthogonality: segments, which are detected as near orthogonal, are made exactly orthogonal.

This is an utility function based on regularize_segments() that is using default parameters.

Template Parameters
InputRangea model of ConstRange whose iterator type is RandomAccessIterator
Parameters
input_rangea const range of input segments for angle regularization
Precondition
input_range.size() >= 2
See also
regularize_segments()
regularize_offsets()

◆ regularize_offsets()

template<typename InputRange >
void CGAL::Shape_regularization::Segments::regularize_offsets ( InputRange &  input_range)

#include <CGAL/Shape_regularization/regularize_segments.h>

regularizes offsets in a set of 2D segments.

Given a set of parallel 2D segments, this function enables to reinforce the collinearity property among these segments that is all parallel segments, which are detected as near collinear, are made exactly collinear.

This is an utility function based on regularize_segments() that is using default parameters.

Template Parameters
InputRangea model of ConstRange whose iterator type is RandomAccessIterator
Parameters
input_rangea const range of input segments for offset regularization
Precondition
input_range.size() >= 2
See also
regularize_segments()
regularize_angles()

◆ regularize_segments()

template<typename InputRange , typename NeighQuery , typename RegType , typename QPSolver , typename NamedParameters = parameters::Default_named_parameters>
void CGAL::Shape_regularization::Segments::regularize_segments ( InputRange &  input_range,
NeighQuery &  neighbor_query,
RegType &  regularization_type,
QPSolver &  quadratic_program,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Shape_regularization/regularize_segments.h>

regularizes a set of 2D segments.

Given a set of unordered 2D segments, this function enables to reinforce three types of regularities among these segments:

  • Parallelism: segments, which are detected as near parallel, are made exactly parallel.
  • Orthogonality: segments, which are detected as near orthogonal, are made exactly orthogonal.
  • Collinearity: parallel segments, which are detected as near collinear, are made exactly collinear.

The user has to provide a NeighborQuery model to access local neighbors of a segment and a RegularizationType model to define the type of regularities that should be addressed. The function is based on the class QP_regularization. Please refer to that class and these concepts for more information.

Template Parameters
InputRangea model of ConstRange whose iterator type is RandomAccessIterator
NeighQuerya model of NeighborQuery
RegTypea model of RegularizationType
QPSolvera model of QuadraticProgramTraits
NamedParametersa sequence of Named Parameters
Parameters
input_rangea const range of input segments for shape regularization
neighbor_queryan instance of NeighQuery that is used internally to access neighbors of a segment; this parameter can be omitted together with the regularization_type parameter, in this case, all types of regularities will be reinforced on the whole input range and the default solver will be used (see below)
regularization_typean instance of RegType that is used internally to obtain bounds and target values required by the regularization; this parameter can be omitted together with the neighbor_query parameter, in this case, all types of regularities will be reinforced on the whole input range and the default solver will be used (see below)
quadratic_programan instance of QPSolver to solve the quadratic programming problem; this parameter can be omitted, the default solver is CGAL::OSQP_quadratic_program_traits
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
Precondition
input_range.size() >= 2
See also
regularize_angles()
regularize_offsets()
Examples:
Shape_regularization/regularize_100_segments_angles.cpp, Shape_regularization/regularize_100_segments_offsets.cpp, Shape_regularization/regularize_real_data_2.cpp, and Shape_regularization/regularize_simple.cpp.

◆ unique_segments()

template<typename InputRange , typename OutIterator , typename NamedParameters = parameters::Default_named_parameters>
OutIterator CGAL::Shape_regularization::Segments::unique_segments ( const InputRange &  input_range,
OutIterator  segments,
const NamedParameters &  np = parameters::default_values() 
)

#include <CGAL/Shape_regularization/regularize_segments.h>

substitutes groups of 2D collinear segments by average segments.

This function 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.

This function does not regularize input segments, but only groups and then simplifies them.

Template Parameters
InputRangea model of ConstRange whose iterator type is RandomAccessIterator
OutIteratora model of OutputIterator that accepts segments of type GeomTraits::Segment_2
NamedParametersa sequence of Named Parameters
Parameters
input_rangea const range of input segments
segmentsan output iterator with the simplified 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
  • maximum allowed orthogonal distance between two parallel segments such that they are considered to be collinear
  • Type: GeomTraits::FT
  • Default: 0.2 unit length
  • indicates whether the order of input segments should be preserved or not
  • Type: boolean
  • Default: false
  • a property map that maps an item from input_range to GeomTraits::Segment_2
  • Type: a model of ReadablePropertyMap whose key type is the value type of the input range and value type is GeomTraits::Segment_2
  • Default: CGAL::Identity_property_map
Returns
an output iterator to the element in the destination range, one past the last segment stored
Precondition
input_range.size() >= 1
maximum_offset >= 0