mirror of https://github.com/CGAL/cgal
add constructor for Polygon_mesh_slicer_3 taking a pre-built AABB_tree of edges
AABB_tree is made a template parameter of Polygon_mesh_slicer_3 also add a test for compilation of this new constructor
This commit is contained in:
parent
af7265f541
commit
ef4df5f47d
|
|
@ -49,18 +49,14 @@ namespace CGAL {
|
||||||
///
|
///
|
||||||
/// Depends on \ref PkgAABB_treeSummary
|
/// Depends on \ref PkgAABB_treeSummary
|
||||||
/// \todo `PolygonMesh` should be a model of `FaceListGraph`
|
/// \todo `PolygonMesh` should be a model of `FaceListGraph`
|
||||||
/// \todo Add a constructor from an AABB-tree (the type is hardcoded given `PolygonMesh`)
|
template<class PolygonMesh,
|
||||||
template<class PolygonMesh, class Kernel>
|
class Kernel,
|
||||||
|
class AABB_tree_ = AABB_tree<
|
||||||
|
AABB_traits<Kernel,
|
||||||
|
AABB_halfedge_graph_segment_primitive<PolygonMesh> > > >
|
||||||
class Polygon_mesh_slicer_3
|
class Polygon_mesh_slicer_3
|
||||||
{
|
{
|
||||||
private:
|
private:
|
||||||
typedef AABB_halfedge_graph_segment_primitive<PolygonMesh> AABB_primitive;
|
|
||||||
typedef AABB_traits<Kernel, AABB_primitive> AABB_traits_;
|
|
||||||
typedef AABB_tree<AABB_traits_> AABB_tree_;
|
|
||||||
|
|
||||||
typedef typename AABB_tree_::Object_and_primitive_id Object_and_primitive_id;
|
|
||||||
typedef typename AABB_tree_::Primitive_id Primitive_id;
|
|
||||||
|
|
||||||
typedef typename Kernel::Plane_3 Plane;
|
typedef typename Kernel::Plane_3 Plane;
|
||||||
typedef typename Kernel::Segment_3 Segment;
|
typedef typename Kernel::Segment_3 Segment;
|
||||||
typedef typename Kernel::Point_3 Point;
|
typedef typename Kernel::Point_3 Point;
|
||||||
|
|
@ -132,7 +128,7 @@ private:
|
||||||
|
|
||||||
// member variables //
|
// member variables //
|
||||||
typename Kernel::Intersect_3 intersect_3_functor;
|
typename Kernel::Intersect_3 intersect_3_functor;
|
||||||
AABB_tree_ tree;
|
const AABB_tree_* tree_ptr;
|
||||||
mutable Node_graph node_graph;
|
mutable Node_graph node_graph;
|
||||||
PolygonMesh& m_pmesh;
|
PolygonMesh& m_pmesh;
|
||||||
|
|
||||||
|
|
@ -244,9 +240,9 @@ private:
|
||||||
{
|
{
|
||||||
node_graph.clear();
|
node_graph.clear();
|
||||||
|
|
||||||
// find out intersecting halfedges (note that tree contains edges only with custom comparator)
|
// find out intersecting halfedges (note that tree_ptr contains edges only with custom comparator)
|
||||||
std::vector<Edge_const_handle> intersected_edges;
|
std::vector<Edge_const_handle> intersected_edges;
|
||||||
tree.all_intersected_primitives(plane, std::back_inserter(intersected_edges));
|
tree_ptr->all_intersected_primitives(plane, std::back_inserter(intersected_edges));
|
||||||
|
|
||||||
// create node graph from segments
|
// create node graph from segments
|
||||||
// each node is associated with multiple edges
|
// each node is associated with multiple edges
|
||||||
|
|
@ -425,12 +421,28 @@ public:
|
||||||
* @param pmesh the polygon mesh to be cut
|
* @param pmesh the polygon mesh to be cut
|
||||||
* @param kernel the kernel
|
* @param kernel the kernel
|
||||||
*/
|
*/
|
||||||
Polygon_mesh_slicer_3(const PolygonMesh& pmesh, const Kernel& kernel = Kernel())
|
Polygon_mesh_slicer_3(const PolygonMesh& pmesh,
|
||||||
|
const Kernel& kernel = Kernel())
|
||||||
: intersect_3_functor(kernel.intersect_3_object()),
|
: intersect_3_functor(kernel.intersect_3_object()),
|
||||||
tree( edges(pmesh).first,
|
|
||||||
edges(pmesh).second,
|
|
||||||
pmesh),
|
|
||||||
m_pmesh(const_cast<PolygonMesh&>(pmesh))
|
m_pmesh(const_cast<PolygonMesh&>(pmesh))
|
||||||
|
{
|
||||||
|
tree_ptr = new AABB_tree_(edges(pmesh).first,
|
||||||
|
edges(pmesh).second,
|
||||||
|
pmesh);
|
||||||
|
}
|
||||||
|
|
||||||
|
/**
|
||||||
|
* Constructor. `pmesh` must be a valid polygon mesh as long as this functor is used.
|
||||||
|
* @param pmesh the polygon mesh to be cut
|
||||||
|
* @param tree a `CGAL::AABB_tree` containing the edges of `pmesh`
|
||||||
|
* @param kernel the kernel
|
||||||
|
*/
|
||||||
|
Polygon_mesh_slicer_3(const PolygonMesh& pmesh,
|
||||||
|
const AABB_tree_& tree,
|
||||||
|
const Kernel& kernel = Kernel())
|
||||||
|
: intersect_3_functor(kernel.intersect_3_object()),
|
||||||
|
tree_ptr(&tree),
|
||||||
|
m_pmesh(const_cast<PolygonMesh&>(pmesh))
|
||||||
{ }
|
{ }
|
||||||
|
|
||||||
/**
|
/**
|
||||||
|
|
@ -445,6 +457,11 @@ public:
|
||||||
CGAL_precondition(!plane.is_degenerate());
|
CGAL_precondition(!plane.is_degenerate());
|
||||||
return intersect_plane(plane, out);
|
return intersect_plane(plane, out);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
~Polygon_mesh_slicer_3()
|
||||||
|
{
|
||||||
|
delete tree_ptr;
|
||||||
|
}
|
||||||
};
|
};
|
||||||
|
|
||||||
}// end of namespace CGAL
|
}// end of namespace CGAL
|
||||||
|
|
|
||||||
|
|
@ -66,6 +66,7 @@ create_single_source_cgal_program( "test_polyhedron_stitching.cpp" )
|
||||||
create_single_source_cgal_program( "orient_polygon_soup_test.cpp" )
|
create_single_source_cgal_program( "orient_polygon_soup_test.cpp" )
|
||||||
create_single_source_cgal_program( "orient_polyhedron_3_test.cpp" )
|
create_single_source_cgal_program( "orient_polyhedron_3_test.cpp" )
|
||||||
create_single_source_cgal_program( "point_inside_polyhedron_boundary_test.cpp" )
|
create_single_source_cgal_program( "point_inside_polyhedron_boundary_test.cpp" )
|
||||||
|
create_single_source_cgal_program( "polygon_mesh_slicer_test.cpp" )
|
||||||
create_single_source_cgal_program( "self_intersection_polyhedron_test.cpp" )
|
create_single_source_cgal_program( "self_intersection_polyhedron_test.cpp" )
|
||||||
create_single_source_cgal_program( "self_intersection_surface_mesh_test.cpp" )
|
create_single_source_cgal_program( "self_intersection_surface_mesh_test.cpp" )
|
||||||
create_single_source_cgal_program( "triangulate_hole_Polyhedron_3_test.cpp" )
|
create_single_source_cgal_program( "triangulate_hole_Polyhedron_3_test.cpp" )
|
||||||
|
|
|
||||||
|
|
@ -0,0 +1,35 @@
|
||||||
|
|
||||||
|
#include <CGAL/Exact_predicates_inexact_constructions_kernel.h>
|
||||||
|
#include <CGAL/Surface_mesh.h>
|
||||||
|
#include <CGAL/boost/graph/graph_traits_Surface_mesh.h>
|
||||||
|
#include <CGAL/AABB_halfedge_graph_segment_primitive.h>
|
||||||
|
|
||||||
|
#include <CGAL/Polygon_mesh_slicer_3.h>
|
||||||
|
#include <CGAL/AABB_tree.h>
|
||||||
|
#include <CGAL/AABB_traits.h>
|
||||||
|
|
||||||
|
#include <fstream>
|
||||||
|
|
||||||
|
typedef CGAL::Exact_predicates_inexact_constructions_kernel K;
|
||||||
|
typedef CGAL::Surface_mesh<K::Point_3> Mesh;
|
||||||
|
|
||||||
|
typedef CGAL::AABB_halfedge_graph_segment_primitive<Mesh> HGSP;
|
||||||
|
typedef CGAL::AABB_traits<K, HGSP> AABB_traits;
|
||||||
|
typedef CGAL::AABB_tree<AABB_traits> AABB_tree;
|
||||||
|
|
||||||
|
int main(int, char** argv)
|
||||||
|
{
|
||||||
|
std::ifstream input("data/U.off");
|
||||||
|
Mesh m;
|
||||||
|
|
||||||
|
if (!input || !(input >> m)){
|
||||||
|
std::cerr << "Error: can not read file.";
|
||||||
|
return 1;
|
||||||
|
}
|
||||||
|
|
||||||
|
AABB_tree tree(edges(m).first, edges(m).second, m);
|
||||||
|
|
||||||
|
CGAL::Polygon_mesh_slicer_3<Mesh, K> slicer(m, tree);
|
||||||
|
|
||||||
|
return 0;
|
||||||
|
}
|
||||||
Loading…
Reference in New Issue