Fixed compilation errors (Linux/g++)

This commit is contained in:
Laurent Saboret 2009-04-27 10:40:34 +00:00
parent f7108fc79f
commit 9c77d89703
6 changed files with 53 additions and 57 deletions

View File

@ -37,10 +37,10 @@ typedef K::Plane_3 Plane;
typedef K::Segment_3 Segment; typedef K::Segment_3 Segment;
typedef K::Triangle_3 Triangle; typedef K::Triangle_3 Triangle;
typedef std::list<typename Segment>::iterator Iterator; typedef std::list<Segment>::iterator Iterator;
typedef CGAL::AABB_segment_primitive<K,Iterator> Primitive; typedef CGAL::AABB_segment_primitive<K,Iterator> Primitive;
typedef CGAL::AABB_traits<typename K, typename Primitive> Traits; typedef CGAL::AABB_traits<K, Primitive> Traits;
typedef CGAL::AABB_tree<typename Traits> Tree; typedef CGAL::AABB_tree<Traits> Tree;
int main(void) int main(void)
{ {
@ -57,11 +57,11 @@ int main(void)
Tree tree(segments.begin(),segments.end()); Tree tree(segments.begin(),segments.end());
Plane plane(a,b,d); Plane plane(a,b,d);
std::cout << tree.number_of_intersections(plane) std::cout << tree.number_of_intersections(plane)
<< " intersections(s) with plane" << std::endl; << " intersections(s) with plane" << std::endl;
Triangle triangle(a,b,c); Triangle triangle(a,b,c);
std::cout << tree.number_of_intersections(triangle) std::cout << tree.number_of_intersections(triangle)
<< " intersections(s) with triangle" << std::endl; << " intersections(s) with triangle" << std::endl;
// TOFIX: following does not compile due to intersection(const Sphere& sphere, // TOFIX: following does not compile due to intersection(const Sphere& sphere,

View File

@ -37,10 +37,10 @@ typedef K::Line_3 Line;
typedef K::Point_3 Point; typedef K::Point_3 Point;
typedef K::Triangle_3 Triangle; typedef K::Triangle_3 Triangle;
typedef std::list<typename Triangle>::iterator Iterator; typedef std::list<Triangle>::iterator Iterator;
typedef CGAL::AABB_triangle_primitive<K,Iterator> Primitive; typedef CGAL::AABB_triangle_primitive<K,Iterator> Primitive;
typedef CGAL::AABB_traits<typename K, typename Primitive> AABB_triangle_traits; typedef CGAL::AABB_traits<K, Primitive> AABB_triangle_traits;
typedef CGAL::AABB_tree<typename AABB_triangle_traits> Tree; typedef CGAL::AABB_tree<AABB_triangle_traits> Tree;
int main(void) int main(void)
{ {
@ -55,16 +55,16 @@ int main(void)
triangles.push_back(Triangle(a,d,c)); triangles.push_back(Triangle(a,d,c));
Tree tree(triangles.begin(),triangles.end()); Tree tree(triangles.begin(),triangles.end());
// counts #intersections // counts #intersections
Ray ray(a,b); Ray ray(a,b);
std::cout << tree.number_of_intersections(ray) std::cout << tree.number_of_intersections(ray)
<< " intersections(s) with ray" << std::endl; << " intersections(s) with ray" << std::endl;
// compute closest point // compute closest point
Point hint(a); Point hint(a);
Point query(2.0, 2.0, 2.0); Point query(2.0, 2.0, 2.0);
Point closest = tree.closest_point(query,hint); Point closest = tree.closest_point(query,hint);
return 0; return 0;
} }

View File

@ -27,6 +27,7 @@
#include <CGAL/Bbox_3.h> #include <CGAL/Bbox_3.h>
#include <CGAL/AABB_intersections.h> #include <CGAL/AABB_intersections.h>
#include <CGAL/Kernel/global_functions.h>
namespace CGAL { namespace CGAL {
@ -50,7 +51,7 @@ public:
/// AABBTraits concept types /// AABBTraits concept types
typedef typename CGAL::Bbox_3 Bounding_box; typedef typename CGAL::Bbox_3 Bounding_box;
typedef typename AABB_primitive Primitive; typedef AABB_primitive Primitive;
typedef typename AABB_primitive::Datum Datum; typedef typename AABB_primitive::Datum Datum;
typedef typename GeomTraits::Sphere_3 Sphere; typedef typename GeomTraits::Sphere_3 Sphere;
@ -136,11 +137,6 @@ public:
bool is_contained(const Sphere& a, const Sphere& b) const; bool is_contained(const Sphere& a, const Sphere& b) const;
private:
/// Private types
typedef typename GeomTraits::FT FT;
typedef typename GeomTraits::Point_3 Point_3;
private: private:
/** /**
* @brief Computes bounding box of one primitive * @brief Computes bounding box of one primitive

View File

@ -47,13 +47,13 @@ namespace CGAL {
typedef typename AABBTraits::Intersection Intersection; typedef typename AABBTraits::Intersection Intersection;
typedef typename AABBTraits::Projection_query Projection_query; typedef typename AABBTraits::Projection_query Projection_query;
private: private:
typedef typename AABB_search_tree<AABBTraits> Search_tree; typedef AABB_search_tree<AABBTraits> Search_tree;
public: public:
/** /**
* @brief Constructor * @brief Constructor
* @param first iterator over first primitive to insert * @param first iterator over first primitive to insert
* @param beyond past-the-end iterator * @param beyond past-the-end iterator
* *
* Builds the datastructure. Type ConstPrimitiveIterator can be any const * Builds the datastructure. Type ConstPrimitiveIterator can be any const
* iterator on a container of Primitive::id_type such that Primitive has * iterator on a container of Primitive::id_type such that Primitive has
@ -348,7 +348,7 @@ namespace CGAL {
{ {
// iterate over primitives to get points on them // iterate over primitives to get points on them
std::list<Projection_query> points; std::list<Projection_query> points;
std::vector<Primitive>::const_iterator it; typename std::vector<Primitive>::const_iterator it;
for(it = m_data.begin(); it != m_data.end(); it++) for(it = m_data.begin(); it != m_data.end(); it++)
{ {
const Primitive& pr = *it; const Primitive& pr = *it;
@ -435,7 +435,7 @@ namespace CGAL {
template<typename Tr> template<typename Tr>
typename AABB_tree<Tr>::Projection typename AABB_tree<Tr>::Projection
AABB_tree<Tr>::closest_point(const Projection_query& query, AABB_tree<Tr>::closest_point(const Projection_query& query,
const Projection& hint) const Projection& hint)
{ {
Projecting_traits traversal_traits(query,hint); Projecting_traits traversal_traits(query,hint);
@ -447,7 +447,7 @@ namespace CGAL {
// first nearest neighbor point to get a hint // first nearest neighbor point to get a hint
template<typename Tr> template<typename Tr>
typename AABB_tree<Tr>::Projection typename AABB_tree<Tr>::Projection
AABB_tree<Tr>::closest_point(const Projection_query& query) AABB_tree<Tr>::closest_point(const Projection_query& query)
{ {
// construct search KD-tree if needed // construct search KD-tree if needed
if(!m_search_tree_constructed) if(!m_search_tree_constructed)

View File

@ -43,14 +43,14 @@ void test_all_query_types(Tree& tree)
{ {
std::cout << "Test all query types" << std::endl; std::cout << "Test all query types" << std::endl;
typedef K::FT FT; typedef typename K::FT FT;
typedef K::Ray_3 Ray; typedef typename K::Ray_3 Ray;
typedef K::Line_3 Line; typedef typename K::Line_3 Line;
typedef K::Point_3 Point; typedef typename K::Point_3 Point;
typedef K::Vector_3 Vector; typedef typename K::Vector_3 Vector;
typedef K::Segment_3 Segment; typedef typename K::Segment_3 Segment;
typedef Tree::Primitive Primitive; typedef typename Tree::Primitive Primitive;
typedef Tree::Intersection Intersection; typedef typename Tree::Intersection Intersection;
Point p((FT)-0.5, (FT)-0.5, (FT)-0.5); Point p((FT)-0.5, (FT)-0.5, (FT)-0.5);
Point q((FT) 0.5, (FT) 0.5, (FT) 0.5); Point q((FT) 0.5, (FT) 0.5, (FT) 0.5);
@ -60,14 +60,14 @@ void test_all_query_types(Tree& tree)
bool success = false; bool success = false;
// do_intersect // do_intersect
success = tree.do_intersect(ray); success = tree.do_intersect(ray);
success = tree.do_intersect(line); success = tree.do_intersect(line);
success = tree.do_intersect(segment); success = tree.do_intersect(segment);
// number_of_intersections // number_of_intersections
tree.number_of_intersections(ray); tree.number_of_intersections(ray);
tree.number_of_intersections(line); tree.number_of_intersections(line);
tree.number_of_intersections(segment); tree.number_of_intersections(segment);
// all_intersected_primitives // all_intersected_primitives
std::list<Primitive> primitives; std::list<Primitive> primitives;
@ -89,14 +89,14 @@ void test_all_query_types(Tree& tree)
} }
template <class Tree, class Polyhedron, class K> template <class Tree, class Polyhedron, class K>
void test_speed(Tree& tree, void test_speed(Tree& tree,
Polyhedron& polyhedron) Polyhedron& polyhedron)
{ {
std::cout << "Test for speed" << std::endl; std::cout << "Test for speed" << std::endl;
typedef K::FT FT; typedef typename K::FT FT;
typedef K::Ray_3 Ray; typedef typename K::Ray_3 Ray;
typedef K::Point_3 Point; typedef typename K::Point_3 Point;
typedef K::Vector_3 Vector; typedef typename K::Vector_3 Vector;
CGAL::Timer timer; CGAL::Timer timer;
unsigned int nb = 0; unsigned int nb = 0;
@ -106,7 +106,7 @@ void test_speed(Tree& tree,
Ray ray(source, vec); Ray ray(source, vec);
while(timer.time() < 1.0) while(timer.time() < 1.0)
{ {
tree.do_intersect(ray); tree.do_intersect(ray);
nb++; nb++;
} }
double speed = (double)nb / timer.time(); double speed = (double)nb / timer.time();
@ -117,10 +117,10 @@ void test_speed(Tree& tree,
template <class K> template <class K>
void test(const char *filename) void test(const char *filename)
{ {
typedef K::FT FT; typedef typename K::FT FT;
typedef K::Ray_3 Ray; typedef typename K::Ray_3 Ray;
typedef K::Point_3 Point; typedef typename K::Point_3 Point;
typedef K::Vector_3 Vector; typedef typename K::Vector_3 Vector;
typedef CGAL::Polyhedron_3<K> Polyhedron; typedef CGAL::Polyhedron_3<K> Polyhedron;
typedef CGAL::AABB_polyhedron_triangle_primitive<K,Polyhedron> Primitive; typedef CGAL::AABB_polyhedron_triangle_primitive<K,Polyhedron> Primitive;
typedef CGAL::AABB_traits<K, Primitive> Traits; typedef CGAL::AABB_traits<K, Primitive> Traits;

View File

@ -33,13 +33,13 @@
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h> #include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
template <class Tree, class Polyhedron, class K> template <class Tree, class Polyhedron, class K>
void test_speed(Tree& tree, void test_speed(Tree& tree,
Polyhedron& polyhedron) Polyhedron& polyhedron)
{ {
typedef K::FT FT; typedef typename K::FT FT;
typedef K::Ray_3 Ray; typedef typename K::Ray_3 Ray;
typedef K::Point_3 Point; typedef typename K::Point_3 Point;
typedef K::Vector_3 Vector; typedef typename K::Vector_3 Vector;
CGAL::Timer timer; CGAL::Timer timer;
unsigned int nb = 0; unsigned int nb = 0;
@ -47,7 +47,7 @@ void test_speed(Tree& tree,
Point query((FT)0.0, (FT)0.0, (FT)0.0); Point query((FT)0.0, (FT)0.0, (FT)0.0);
while(timer.time() < 1.0) while(timer.time() < 1.0)
{ {
Point closest = tree.closest_point(query); Point closest = tree.closest_point(query);
nb++; nb++;
} }
double speed = (double)nb / timer.time(); double speed = (double)nb / timer.time();
@ -58,10 +58,10 @@ void test_speed(Tree& tree,
template <class K> template <class K>
void test(const char *filename) void test(const char *filename)
{ {
typedef K::FT FT; typedef typename K::FT FT;
typedef K::Ray_3 Ray; typedef typename K::Ray_3 Ray;
typedef K::Point_3 Point; typedef typename K::Point_3 Point;
typedef K::Vector_3 Vector; typedef typename K::Vector_3 Vector;
typedef CGAL::Polyhedron_3<K> Polyhedron; typedef CGAL::Polyhedron_3<K> Polyhedron;
typedef CGAL::AABB_polyhedron_triangle_primitive<K,Polyhedron> Primitive; typedef CGAL::AABB_polyhedron_triangle_primitive<K,Polyhedron> Primitive;
typedef CGAL::AABB_traits<K, Primitive> Traits; typedef CGAL::AABB_traits<K, Primitive> Traits;