fix the concept ConvexHullTraits_3 in the case of coplanar points

This commit is contained in:
Sébastien Loriot 2014-06-10 22:42:50 +02:00
parent 4ed2497f2e
commit 7b05330a75
4 changed files with 71 additions and 69 deletions

View File

@ -77,7 +77,11 @@ typedef R::Construct_centroid_3 Construct_centroid_3;
/*!
*/
typedef unspecified_type Construct_orthogonal_vector_3;
typedef R::Construct_vector_3 Construct_vector_3;
/*!
*/
typedef R::Orientation_3 Orientation_3;
/*!

View File

@ -128,8 +128,25 @@ Same as above but in the `xz`-plane
*/
typedef unspecified_type Traits_xz_3;
/*!
Function object type that provides
`Vector_3 operator()(Point_3 p1, Point_3 p2)`, which
constructs and returns the vector `p1p2`, and
`Vector_3 operator()(int x,int y,int z)` which constructs and returns
a vector with integer Cartesian coordinates.
*/
typedef unspecified_type Construct_vector_3;
/// @}
/*!
Predicate object type that
provides `CGAL::Orientation operator()(Vector_3 u, Vector_3 v, Vector_3 w)`
which returns `CGAL::NEGATIVE` if u, v and w are negatively oriented,
`CGAL::POSITIVE` if u, v and w are positively oriented and
`CGAL::COPLANAR` if u, v and w are coplanar.
*/
typedef unspecified_type Orientation_3;
/// @}
/// \name Creation
/// Only a copy constructor is required.

View File

@ -130,29 +130,6 @@ public:
}
};
template <class T>
class Max_coordinate_3
{
public:
int operator()(const T& v)
{
if (CGAL_NTS abs(v.x()) >= CGAL_NTS abs(v.y()))
{
if (CGAL_NTS abs(v.x()) >= CGAL_NTS abs(v.z())) return 0;
return 2;
}
else
{
if (CGAL_NTS abs(v.y()) >= CGAL_NTS abs(v.z())) return 1;
return 2;
}
}
};
template <typename GT>
struct GT3_for_CH3 {
typedef typename GT::Point_3 Point_2;
@ -198,13 +175,13 @@ class Convex_hull_traits_3
typedef Point_triple_has_on_positive_side_3<Self> Has_on_positive_side_3;
typedef Point_triple_less_signed_distance_to_plane_3<Self, R>
Less_signed_distance_to_plane_3;
Less_signed_distance_to_plane_3;
// required for degenerate case of all points coplanar
typedef CGAL::Max_coordinate_3<Vector_3> Max_coordinate_3;
typedef CGAL::Projection_traits_xy_3<R> Traits_xy_3;
typedef CGAL::Projection_traits_yz_3<R> Traits_yz_3;
typedef CGAL::Projection_traits_xz_3<R> Traits_xz_3;
typedef CGAL::Projection_traits_xy_3<R> Traits_xy_3;
typedef CGAL::Projection_traits_yz_3<R> Traits_yz_3;
typedef CGAL::Projection_traits_xz_3<R> Traits_xz_3;
typedef typename R::Construct_vector_3 Construct_vector_3;
// for postcondition checking
typedef typename R::Ray_3 Ray_3;
@ -272,9 +249,14 @@ class Convex_hull_traits_3
less_signed_distance_to_plane_3_object() const
{ return Less_signed_distance_to_plane_3(); }
Max_coordinate_3
max_coordinate_3_object() const
{ return Max_coordinate_3(); }
Orientation_3
orientation_3_object() const
{ return Orientation_3(); }
Construct_vector_3
construct_vector_3_object() const
{ return Construct_vector_3(); }
};
} // namespace CGAL

View File

@ -309,14 +309,15 @@ struct Projection_traits<T,true>{
} } //end of namespace internal::Convex_hull_3
template <class InputIterator, class Plane_3, class Polyhedron_3, class Traits>
template <class InputIterator, class Point_3, class Polyhedron_3, class Traits>
void coplanar_3_hull(InputIterator first, InputIterator beyond,
Plane_3 plane, Polyhedron_3& P, const Traits& traits)
const Point_3& p1, const Point_3& p2, const Point_3& p3,
Polyhedron_3& P, const Traits& traits)
{
typedef typename Traits::Point_3 Point_3;
typedef typename Traits::Vector_3 Vector_3;
typedef Max_coordinate_3<Vector_3> Max_coordinate_3;
typedef Polyhedron_3 Polyhedron;
typedef typename Traits::Construct_vector_3 Construct_vector_3;
typedef typename Traits::Orientation_3 Orientation_3;
typedef typename internal::Convex_hull_3::Projection_traits<Traits> PTraits;
typedef typename PTraits::Traits_xy_3 Traits_xy_3;
typedef typename PTraits::Traits_yz_3 Traits_yz_3;
@ -324,38 +325,36 @@ void coplanar_3_hull(InputIterator first, InputIterator beyond,
std::list<Point_3> CH_2;
typedef typename std::list<Point_3>::iterator CH_2_iterator;
typedef typename Traits::Construct_orthogonal_vector_3
Construct_normal_vec;
Max_coordinate_3 max_coordinate;
Construct_vector_3 vector_3 = traits.construct_vector_3_object();
Orientation_3 orientation = traits.orientation_3_object();
Vector_3 v1 = vector_3(p1,p2);
Vector_3 v2 = vector_3(p1,p3);
Vector_3 vx = vector_3(1,0,0);
Construct_normal_vec c_normal =
traits.construct_orthogonal_vector_3_object();
Vector_3 normal = c_normal(plane);
int max_coord = max_coordinate(normal);
switch (max_coord)
{
case 0:
{
convex_hull_points_2(first, beyond,
std::back_inserter(CH_2), Traits_yz_3());
break;
}
case 1:
{
convex_hull_points_2(first, beyond,
std::back_inserter(CH_2), Traits_xz_3());
break;
}
case 2:
{
convex_hull_points_2(first, beyond,
std::back_inserter(CH_2), Traits_xy_3());
break;
}
default:
break;
if ( orientation(v1, v2, vx) != COPLANAR )
convex_hull_points_2( first, beyond,
std::back_inserter(CH_2),
Traits_yz_3() );
else{
Vector_3 vy = vector_3(0,1,0);
if ( orientation(v1,v2,vy) != COPLANAR )
convex_hull_points_2( first, beyond,
std::back_inserter(CH_2),
Traits_xz_3() );
else{
CGAL_assertion_code( Vector_3 vz = vector_3(0,0,1); )
CGAL_assertion( orientation(v1,v2,vz) != COPLANAR );
convex_hull_points_2( first, beyond,
std::back_inserter(CH_2),
Traits_xy_3() );
}
}
typedef typename Polyhedron::Halfedge_data_structure HDS;
typedef typename Polyhedron_3::Halfedge_data_structure HDS;
Build_coplanar_poly<HDS,CH_2_iterator> poly(CH_2.begin(),CH_2.end());
P.delegate(poly);
@ -745,7 +744,7 @@ ch_quickhull_polyhedron_3(std::list<typename Traits::Point_3>& points,
// if the maximum distance point is on the plane then all are coplanar
if (coplanar(*point1_it, *point2_it, *point3_it, *max_it)) {
coplanar_3_hull(points.begin(), points.end(), plane, P, traits);
coplanar_3_hull(points.begin(), points.end(), *point1_it, *point2_it, *point3_it, P, traits);
} else {
Tds tds;
Vertex_handle v0 = tds.create_vertex(); v0->set_point(*point1_it);