Add Default to Convex_hull_graph_traits_3; BGLize some code

This commit is contained in:
Andreas Fabri 2016-11-22 12:27:18 +01:00
parent d968e123ea
commit e26e60dd6c
7 changed files with 35 additions and 23 deletions

View File

@ -29,6 +29,7 @@
#include <list>
#include <CGAL/Filtered_predicate.h>
#include <CGAL/Cartesian_converter.h>
#include <CGAL/Default.h>
namespace CGAL {
template < class R_ >
@ -185,7 +186,7 @@ struct Convex_hull_traits_base_3<R_, Tag_true>{
};
template <class R_, class Polyhedron = Polyhedron_3<R_>, class Has_filtered_predicates_tag = Tag_false>
template <class R_, class Polyhedron = Default, class Has_filtered_predicates_tag = Tag_false>
class Convex_hull_traits_3 :
public Convex_hull_traits_base_3<R_, Has_filtered_predicates_tag>
{
@ -198,7 +199,7 @@ class Convex_hull_traits_3 :
typedef Point_triple<R> Plane_3;
typedef typename R::Vector_3 Vector_3;
typedef Polyhedron Polyhedron_3;
typedef typename Default::Get<Polyhedron, Polyhedron_3<R> >::type Polyhedron_3;
typedef typename R::Construct_segment_3 Construct_segment_3;
typedef typename R::Construct_ray_3 Construct_ray_3;

View File

@ -71,7 +71,7 @@ struct Default_traits_for_Chull_3{
//FT is a floating point type and Kernel is a filtered kernel
template <class Point_3>
struct Default_traits_for_Chull_3<Point_3,boost::true_type,Tag_true>{
typedef Convex_hull_traits_3< typename Kernel_traits<Point_3>::Kernel, Tag_true > type;
typedef Convex_hull_traits_3< typename Kernel_traits<Point_3>::Kernel, Default, Tag_true > type;
};
template <class Traits>
@ -79,9 +79,9 @@ struct Default_polyhedron_for_Chull_3{
typedef CGAL::Polyhedron_3<Traits> type;
};
template <class K,class Tag>
struct Default_polyhedron_for_Chull_3<Convex_hull_traits_3<K, Tag> >{
typedef typename Convex_hull_traits_3<K, Tag>::Polyhedron_3 type;
template <class K, class P, class Tag>
struct Default_polyhedron_for_Chull_3<Convex_hull_traits_3<K, P, Tag> >{
typedef typename Convex_hull_traits_3<K, P, Tag>::Polyhedron_3 type;
};
//utility class to select the right version of internal predicate Is_on_positive_side_of_plane_3
@ -130,11 +130,11 @@ public:
//The main operator() first tries the static version of the predicate, then uses
//interval arithmetic (the protector must be created before using this predicate)
//and in case of failure, exact arithmetic is used.
template <class Kernel>
class Is_on_positive_side_of_plane_3<Convex_hull_traits_3<Kernel, Tag_true>,Tag_true>{
typedef Simple_cartesian<CGAL::internal::Exact_field_selector<double>::Type> PK;
template <class Kernel, class P>
class Is_on_positive_side_of_plane_3<Convex_hull_traits_3<Kernel, P, Tag_true>,Tag_true>{
typedef Simple_cartesian<CGAL::internal::Exact_field_selector<double>::Type> PK;
typedef Simple_cartesian<Interval_nt_advanced > CK;
typedef Convex_hull_traits_3<Kernel, Tag_true> Traits;
typedef Convex_hull_traits_3<Kernel, P, Tag_true> Traits;
typedef typename Traits::Point_3 Point_3;
Cartesian_converter<Kernel,CK> to_CK;
@ -616,6 +616,7 @@ ch_quickhull_polyhedron_3(std::list<typename Traits::Point_3>& points,
typedef Triangulation_data_structure_2<
Triangulation_vertex_base_with_info_2<int, GT3_for_CH3<Traits> >,
Convex_hull_face_base_2<int, Traits> > Tds;
typedef typename Tds::Vertex_handle Vertex_handle;
typedef typename Tds::Face_handle Face_handle;
@ -685,6 +686,7 @@ ch_quickhull_polyhedron_3(std::list<typename Traits::Point_3>& points,
} } //namespace internal::Convex_hull_3
template <class InputIterator, class Traits>
void
convex_hull_3(InputIterator first, InputIterator beyond,
@ -768,7 +770,8 @@ convex_hull_3(InputIterator first, InputIterator beyond,
}
// result will be a polyhedron
typename internal::Convex_hull_3::Default_polyhedron_for_Chull_3<Traits>::type P;
typedef typename internal::Convex_hull_3::Default_polyhedron_for_Chull_3<Traits>::type Polyhedron;
Polyhedron P;
P3_iterator minx, maxx, miny, it;
minx = maxx = miny = it = points.begin();
@ -783,14 +786,18 @@ convex_hull_3(InputIterator first, InputIterator beyond,
} else {
internal::Convex_hull_3::ch_quickhull_polyhedron_3(points, point1_it, point2_it, point3_it, P, traits);
}
CGAL_assertion(P.size_of_vertices()>=3);
if (boost::next(P.vertices_begin(),3) == P.vertices_end()){
CGAL_assertion(num_vertices(P)>=3);
boost::graph_traits<Polyhedron>::halfedge_iterator b,e;
boost::tie(b,e) = halfedges(P);
if (boost::next(b,3) == e){
typename boost::property_map<Polyhedron, vertex_point_t>::type vpmap = get(CGAL::vertex_point, P);
typedef typename Traits::Triangle_3 Triangle_3;
typename Traits::Construct_triangle_3 construct_triangle =
traits.construct_triangle_3_object();
Triangle_3 tri = construct_triangle(P.halfedges_begin()->vertex()->point(),
P.halfedges_begin()->next()->vertex()->point(),
P.halfedges_begin()->opposite()->vertex()->point());
Triangle_3 tri = construct_triangle(get(vpmap, target(*b,P)),
get(vpmap, target(next(*b,P),P)),
get(vpmap, target(opposite(*b,P),P)));
ch_object = make_object(tri);
}
else

View File

@ -87,7 +87,10 @@ bool is_strongly_convex_3(const Polyhedron& P, const Traits& traits)
typename boost::property_map<Polyhedron, vertex_point_t>::const_type vpmap = get(CGAL::vertex_point, P);
if (P.vertices_begin() == P.vertices_end()) return false;
vertex_iterator v_it, v_it_e;
boost::tie(v_it, v_it_e) = vertices(P);
if (v_it == v_it_e) return false;
BOOST_FOREACH(face_descriptor fd , faces(P))
if (!is_locally_convex(P, vpmap, fd, traits))
@ -98,9 +101,6 @@ bool is_strongly_convex_3(const Polyhedron& P, const Traits& traits)
typename Traits::Coplanar_3 coplanar = traits.coplanar_3_object();
vertex_iterator v_it, v_it_e;
boost::tie(v_it, v_it_e) = vertices(P);
face_iterator f_it, f_it_e;
boost::tie(f_it, f_it_e) = faces(P);
Point_3 p;

View File

@ -25,7 +25,7 @@ int main()
CGAL_static_assertion( (boost::is_same<SCD,Default_traits_for_Chull_3<SCD::Point_3>::type>::value) );
CGAL_static_assertion( (boost::is_same<SCR,Default_traits_for_Chull_3<SCR::Point_3>::type>::value) );
CGAL_static_assertion( (boost::is_same<EPEC,Default_traits_for_Chull_3<EPEC::Point_3>::type>::value) );
CGAL_static_assertion( (boost::is_same<CGAL::Convex_hull_traits_3<EPIC, CGAL::Tag_true>,Default_traits_for_Chull_3<EPIC::Point_3>::type>::value) );
CGAL_static_assertion( (boost::is_same<Is_on_positive_side_of_plane_3<CGAL::Convex_hull_traits_3<EPIC, CGAL::Tag_true> >::Protector,CGAL::Protect_FPU_rounding<true> >::value) );
CGAL_static_assertion( (boost::is_same<CGAL::Convex_hull_traits_3<EPIC, CGAL::Default, CGAL::Tag_true>,Default_traits_for_Chull_3<EPIC::Point_3>::type>::value) );
CGAL_static_assertion( (boost::is_same<Is_on_positive_side_of_plane_3<CGAL::Convex_hull_traits_3<EPIC, CGAL::Default, CGAL::Tag_true> >::Protector,CGAL::Protect_FPU_rounding<true> >::value) );
return 0;
}

View File

@ -151,6 +151,9 @@ void test_collinear()
}
#include <CGAL/Triangulation_face_base_with_info_2.h>
int main()
{
std::vector<Point_3> points;

View File

@ -31,6 +31,7 @@
namespace CGAL {
//mechanism to abuse Handle_hash_function which is the default
//template parameter of Unique_hash_map
namespace internal{

View File

@ -32,7 +32,7 @@ namespace CGAL {
template<class Triangulation_3,class FG>
typename boost::graph_trait<FG>::vertex_descriptor
typename boost::graph_traits<FG>::vertex_descriptor
star_to_face_graph(const Triangulation_3& t,
typename Triangulation_3::Vertex_handle vh,
FG& fg,