CGAL 6.1 - 3D Convex Hulls
Loading...
Searching...
No Matches
Intersection Test Functions

Functions

template<class PointRange , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
bool CGAL::Convex_hull_3::do_intersect (const PointRange &r1, const PointRange &r2, const NamedParameters_1 &np1=parameters::default_values(), const NamedParameters_2 &np2=parameters::default_values())
 checks if the convex hulls of the two point sets intersect or not.
 
template<class AdjacencyGraph , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
bool CGAL::Convex_hull_3::do_intersect (const AdjacencyGraph &g1, const AdjacencyGraph &g2, const NamedParameters_1 &np1=parameters::default_values(), const NamedParameters_2 &np2=parameters::default_values())
 checks if the two convex graphs intersect or not.
 
template<class PolygonMesh , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
bool CGAL::Convex_hull_3::do_intersect (const Convex_hull_hierarchy< PolygonMesh > &ch1, const Convex_hull_hierarchy< PolygonMesh > &ch2, const NamedParameters_1 &np1=parameters::default_values(), const NamedParameters_2 &np2=parameters::default_values())
 checks if the convex hulls intersect or not.
 
template<class PointRange , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
FT CGAL::Convex_hull_3::separation_distance (const PointRange &r1, const PointRange &r2, const NamedParameters_1 &np1=parameters::default_values(), const NamedParameters_2 &np2=parameters::default_values())
 provides a lower bound on the squared distance between the convex hulls of the two point sets.
 
template<class AdjacencyGraph , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
FT CGAL::Convex_hull_3::separation_distance (const AdjacencyGraph &g1, const AdjacencyGraph &g2, const NamedParameters_1 &np1=parameters::default_values(), const NamedParameters_2 &np2=parameters::default_values())
 provides a lower bound on the squared distance between the two convex graphs.
 
template<class PolygonMesh , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
FT CGAL::Convex_hull_3::separation_distance (const Convex_hull_hierarchy< PolygonMesh > &ch1, const Convex_hull_hierarchy< PolygonMesh > &ch2, const NamedParameters_1 &np1=parameters::default_values(), const NamedParameters_2 &np2=parameters::default_values())
 provides a lower bound on the squared distance between the two convex hulls.
 

Function Documentation

◆ do_intersect() [1/3]

template<class AdjacencyGraph , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
bool CGAL::Convex_hull_3::do_intersect ( const AdjacencyGraph g1,
const AdjacencyGraph g2,
const NamedParameters_1 &  np1 = parameters::default_values(),
const NamedParameters_2 &  np2 = parameters::default_values() 
)

#include <CGAL/Convex_hull_3/predicates.h>

checks if the two convex graphs intersect or not.

Template Parameters
AdjacencyGraphis a model of AdjacencyGraph.
NamedParameters_1a sequence of Named Parameters
NamedParameters_2a sequence of Named Parameters
Parameters
g1the first convex graph
g2the second convex graph
np1an optional sequence of Named Parameters among the ones listed below
np2an optional sequence of Named Parameters among the ones listed below
Precondition
The input graph must represent a convex object to guarantee a correct answer.
Optional Named Parameters
  • a property map associating points to the vertices of g1 (g2)
  • Type: a model of ReadablePropertyMap whose value types are the same for np1 and np2
  • Default: boost::get(CGAL::vertex_point, g)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in AdjacencyGraph.
  • An instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: np1 only
  • if not 0 (no limit), the maximum number of iterations performed by the algorithm. If this value is not 0, then an intersection might be reported even if the convex hulls do not intersect. However, if the convex hulls are reported not to intersect, this is guaranteed.
  • Type: a positive integer convertible to std::size_t
  • Extra: np1 only
  • Default: 0

◆ do_intersect() [2/3]

template<class PolygonMesh , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
bool CGAL::Convex_hull_3::do_intersect ( const Convex_hull_hierarchy< PolygonMesh > &  ch1,
const Convex_hull_hierarchy< PolygonMesh > &  ch2,
const NamedParameters_1 &  np1 = parameters::default_values(),
const NamedParameters_2 &  np2 = parameters::default_values() 
)

#include <CGAL/Convex_hull_3/predicates.h>

checks if the convex hulls intersect or not.

Template Parameters
PolygonMeshis a model of MutableFaceGraph, more details in CGAL::Convex_hull_hierarchy
NamedParameters_1a sequence of Named Parameters
NamedParameters_2a sequence of Named Parameters
Parameters
ch1the first convex hull
ch2the second convex hull
np1an optional sequence of Named Parameters among the ones listed below
np2an optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of ch1 (ch2)
  • Type: a model of ReadablePropertyMap whose value types are the same for np1 and np2
  • Default: boost::get(CGAL::vertex_point, g)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in IncidentGraph.
  • An instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: np1 only
  • if not 0 (no limit), indicates the maximum number of iterations performed by the algorithm. If this value is not 0, then an intersection might be reported even if the convex hulls does not intersect. However, if the convex hulls are reported not to intersect, this is guaranteed.
  • Type: a positive integer convertible to std::size_t
  • Extra: np1 only
  • Default: 0

◆ do_intersect() [3/3]

template<class PointRange , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
bool CGAL::Convex_hull_3::do_intersect ( const PointRange &  r1,
const PointRange &  r2,
const NamedParameters_1 &  np1 = parameters::default_values(),
const NamedParameters_2 &  np2 = parameters::default_values() 
)

#include <CGAL/Convex_hull_3/predicates.h>

checks if the convex hulls of the two point sets intersect or not.

Template Parameters
PointRangeis a model of ConstRange. The value type of its iterator is the key type of the named parameter point_map.
NamedParameters_1a sequence of Named Parameters
NamedParameters_2a sequence of Named Parameters
Parameters
r1range points of the first convex considered in the do-intersect test
r2range points of the second convex considered in the do-intersect test
np1an optional sequence of Named Parameters among the ones listed below
np2an optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • An instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: np1 only
  • if not 0 (no limit), the maximum number of iterations performed by the algorithm. If this value is not 0, then an intersection might be reported even if the convex hulls do not intersect. However, if the convex hulls are reported not to intersect, this is guaranteed.
  • Type: a positive integer convertible to std::size_t
  • Extra: np1 only
  • Default: 0
Examples
Convex_hull_3/do_intersect_ch3.cpp.

◆ separation_distance() [1/3]

template<class AdjacencyGraph , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
FT CGAL::Convex_hull_3::separation_distance ( const AdjacencyGraph g1,
const AdjacencyGraph g2,
const NamedParameters_1 &  np1 = parameters::default_values(),
const NamedParameters_2 &  np2 = parameters::default_values() 
)

#include <CGAL/Convex_hull_3/predicates.h>

provides a lower bound on the squared distance between the two convex graphs.

Template Parameters
AdjacencyGraphis a model of AdjacencyGraph.
NamedParameters_1a sequence of Named Parameters
NamedParameters_2a sequence of Named Parameters
Parameters
g1the first convex graph
g2the second convex graph
np1an optional sequence of Named Parameters among the ones listed below
np2an optional sequence of Named Parameters among the ones listed below
Warning
The input graph must represent a convex object to guarantee a correct answer.
Optional Named Parameters
  • a property map associating points to the vertices of g1 (g2)
  • Type: a model of ReadablePropertyMap whose value types are the same for np1 and np2
  • Default: boost::get(CGAL::vertex_point, g)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in AdjacencyGraph.
  • An instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: np1 only
  • if not 0 (no limit), indicates the maximum number of iterations performed by the algorithm. if this value is not 0, then the return value can be zero even if the convex hulls does not intersect. However, the value reported remains a lower bound of the distance between the convex.
  • Type: a positive integer convertible to std::size_t
  • Extra: np1 only
  • Default: 0

◆ separation_distance() [2/3]

template<class PolygonMesh , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
FT CGAL::Convex_hull_3::separation_distance ( const Convex_hull_hierarchy< PolygonMesh > &  ch1,
const Convex_hull_hierarchy< PolygonMesh > &  ch2,
const NamedParameters_1 &  np1 = parameters::default_values(),
const NamedParameters_2 &  np2 = parameters::default_values() 
)

#include <CGAL/Convex_hull_3/predicates.h>

provides a lower bound on the squared distance between the two convex hulls.

Template Parameters
PolygonMeshis a model of MutableFaceGraph, more details in CGAL::Convex_hull_hierarchy
NamedParameters_1a sequence of Named Parameters
NamedParameters_2a sequence of Named Parameters
Parameters
ch1the first convex hull
ch2the second convex hull
np1an optional sequence of Named Parameters among the ones listed below
np2an optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • a property map associating points to the vertices of ch1 (ch2)
  • Type: a model of ReadablePropertyMap whose value types are the same for np1 and np2
  • Default: boost::get(CGAL::vertex_point, g)
  • Extra: If this parameter is omitted, an internal property map for CGAL::vertex_point_t must be available in IncidenceGraph.
  • An instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: np1 only
  • if not 0 (no limit), indicates the maximum number of iterations performed by the algorithm. If this value is not 0, then the return value can be zero even if the convex hulls does not intersect. However, the value reported remains a lower bound of the distance between the convex.
  • Type: a positive integer convertible to std::size_t
  • Extra: np1 only
  • Default: 0

◆ separation_distance() [3/3]

template<class PointRange , class NamedParameters_1 = parameters::Default_named_parameters, class NamedParameters_2 = parameters::Default_named_parameters>
FT CGAL::Convex_hull_3::separation_distance ( const PointRange &  r1,
const PointRange &  r2,
const NamedParameters_1 &  np1 = parameters::default_values(),
const NamedParameters_2 &  np2 = parameters::default_values() 
)

#include <CGAL/Convex_hull_3/predicates.h>

provides a lower bound on the squared distance between the convex hulls of the two point sets.

Template Parameters
PointRangeis a model of ConstRange. The value type of its iterator is the key type of the named parameter point_map.
NamedParameters_1a sequence of Named Parameters
NamedParameters_2a sequence of Named Parameters
Parameters
r1first point range
r2second point range
np1an optional sequence of Named Parameters among the ones listed below
np2an optional sequence of Named Parameters among the ones listed below
Optional Named Parameters
  • An instance of a geometric traits class
  • Type: a class model of Kernel
  • Default: a CGAL kernel deduced from the point type, using CGAL::Kernel_traits
  • Extra: np1 only
  • if not 0 (no limit), indicates the maximum number of iterations performed by the algorithm. If this value is not 0, then an intersection might be reported even if the convex hulls does not intersect. However, if the convex hulls are reported not to intersect, this is guaranteed.
  • Type: a positive integer convertible to std::size_t
  • Extra: np1 only
  • Default: 0