mirror of https://github.com/CGAL/cgal
AABB tree now linked to oracle
This commit is contained in:
parent
dfab7987ac
commit
7032b35789
|
|
@ -38,7 +38,7 @@ void MainWindow::on_actionRemeshing_triggered()
|
|||
|
||||
// AABB tree
|
||||
std::cout << "Build AABB tree...";
|
||||
typedef AABB_tree<Kernel,Polyhedron::Facet_handle,Polyhedron> Tree;
|
||||
typedef CGAL::AABB_tree<GT,Polyhedron::Facet_handle,Polyhedron> Tree;
|
||||
Tree tree;
|
||||
tree.build_faces(*pMesh);
|
||||
std::cout << "done" << std::endl;
|
||||
|
|
@ -46,11 +46,13 @@ void MainWindow::on_actionRemeshing_triggered()
|
|||
// input surface
|
||||
typedef CGAL::AABB_polyhedral_oracle<Polyhedron,GT> Input_surface;
|
||||
|
||||
// instantiate surface (linked to the AABB tree)
|
||||
Input_surface input(&tree);
|
||||
|
||||
// remesh
|
||||
Input_surface input;
|
||||
//std::cout << "Remesh...";
|
||||
std::cout << "Remesh...";
|
||||
CGAL::make_surface_mesh(c2t3, input, facets_criteria, CGAL::Manifold_tag());
|
||||
//std::cout << "done (" << tr.number_of_vertices() << " vertices)" << std::endl;
|
||||
std::cout << "done (" << tr.number_of_vertices() << " vertices)" << std::endl;
|
||||
|
||||
// add remesh as new polyhedron
|
||||
Polyhedron *pRemesh = new Polyhedron;
|
||||
|
|
|
|||
|
|
@ -13,6 +13,7 @@
|
|||
#include "Plucker_ray_3_Bbox_3_do_intersect.h"
|
||||
#include "Sphere_3_Bbox_do_intersect.h"
|
||||
|
||||
CGAL_BEGIN_NAMESPACE
|
||||
|
||||
template <class Kernel, class Input, class PSC>
|
||||
class AABB_node
|
||||
|
|
@ -1512,5 +1513,7 @@ public:
|
|||
}
|
||||
};
|
||||
|
||||
CGAL_END_NAMESPACE
|
||||
|
||||
#endif // _AABB_NODE_
|
||||
|
||||
|
|
|
|||
|
|
@ -5,6 +5,8 @@
|
|||
#include <utility>
|
||||
#include <CGAL/iterator.h>
|
||||
|
||||
#include "AABB_tree.h"
|
||||
|
||||
namespace CGAL {
|
||||
|
||||
template <class Polyhedron, class Kernel>
|
||||
|
|
@ -17,16 +19,33 @@ public:
|
|||
typedef typename Kernel::Line_3 Line_3;
|
||||
typedef typename Kernel::Point_3 Point_3;
|
||||
typedef typename Kernel::Segment_3 Segment_3;
|
||||
|
||||
typedef typename AABB_polyhedral_oracle<Polyhedron,Kernel> Self;
|
||||
typedef typename Self Surface_mesher_traits_3;
|
||||
typedef typename Point_3 Intersection_point;
|
||||
|
||||
public:
|
||||
// AABB tree
|
||||
typedef AABB_tree<Kernel,typename Polyhedron::Facet_handle,Polyhedron> Tree;
|
||||
Tree *m_pTree;
|
||||
|
||||
public:
|
||||
Tree* tree() { return m_pTree; }
|
||||
|
||||
public:
|
||||
// Surface constructor
|
||||
AABB_polyhedral_oracle()
|
||||
{
|
||||
m_pTree = NULL;
|
||||
}
|
||||
AABB_polyhedral_oracle(Tree *pTree)
|
||||
{
|
||||
m_pTree = pTree;
|
||||
}
|
||||
AABB_polyhedral_oracle(const AABB_polyhedral_oracle& oracle)
|
||||
{
|
||||
m_pTree = oracle.tree();
|
||||
}
|
||||
|
||||
|
||||
class Intersect_3;
|
||||
|
||||
|
|
|
|||
|
|
@ -8,6 +8,8 @@
|
|||
#include "AABB_node.h"
|
||||
#include "knn.h"
|
||||
|
||||
CGAL_BEGIN_NAMESPACE
|
||||
|
||||
template <class Kernel, class Input, class PSC>
|
||||
class AABB_tree
|
||||
{
|
||||
|
|
@ -481,4 +483,6 @@ public:
|
|||
}
|
||||
};
|
||||
|
||||
CGAL_END_NAMESPACE
|
||||
|
||||
#endif // _AABB_TREE_
|
||||
|
|
|
|||
Loading…
Reference in New Issue