Subject: CGAL users discussion list
List archive
- From: Peter Hachenberger <>
- To:
- Subject: Re: [cgal-discuss][Nef_polehedron]convert_to_Polyhedron ( Polyhedron& P)
- Date: Thu, 17 Jan 2008 14:12:50 +0100
Hi Max,
I got this question frequently off late. For this reason I implemented a
new function. Add the attached file shell_to_nef_3.h to your
include/CGAL/Nef_3 directory and add the following code to
Nef_polyhedron_3.h:
Nef_polyhedron_3(const Nef_polyhedron& N,
SFace_const_iterator sf)
{
SNC_structure rsnc;
*this = Nef_polyhedron_3(rsnc, new SNC_point_locator_default, false);
initialize_infibox_vertices(EMPTY);
shell_to_nef_3(N, sf, snc());
build_external_structure();
simplify();
CGAL::Mark_bounded_volumes<Nef_polyhedron_3> mbv(true);
delegate(mbv);
set_snc(snc());
}
There is another file attached, which demonstrates the usage. This
program compares two methods of conversion, an old and simple one with
the new one. If you don't use CGAL 3.3.1, you won't have the code for
the old implementation and need to comment/ignore the first part of the
program. I originally wrote the code to work for closed inner shells,
but maybe it also works non-closed shells and outer shells. Please tell
me how it works.
Peter
On Thu, 2008-01-17 at 19:18 +0800, Max wrote:
> Hi all,
>
> I know if a Nef_polyhedron is simple (or, in other words,
> a manifold), it could be converted to a single polyhedron
> with the method:
>
> Nef_polyhedron::convert_to_Polyhedron ( Polyhedron& P)
>
> My question is, if N is not a single manifold, but still finite,
> is there a way in CGAL to convert it into a group of polyhedra?
>
> Thanks for any help.
>
> B/Rgds
> Max
>
// Copyright (c) 1997-2002 Max-Planck-Institute Saarbruecken (Germany). // All rights reserved. // // This file is part of CGAL (www.cgal.org); you may redistribute it under // the terms of the Q Public License version 1.0. // See the file LICENSE.QPL distributed with CGAL. // // Licensees holding a valid commercial license may use this file in // accordance with the commercial license agreement provided with the software. // // This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE // WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE. // // $URL: svn+ssh:///svn/cgal/trunk/Nef_3/include/CGAL/Nef_3/polyhedron_3_to_nef_3.h $ // $Id: polyhedron_3_to_nef_3.h 39613 2007-07-31 11:41:47Z hachenb $ // // // Author(s) : Michael Seel <> // Miguel Granados <> // Susan Hert <> // Lutz Kettner <> #ifndef CGAL_NEF_SHELL_TO_NEF_3_H #define CGAL_NEF_SHELL_TO_NEF_3_H #include <CGAL/Circulator_project.h> #include <CGAL/normal_vector_newell_3.h> #include <CGAL/Nef_S2/SM_point_locator.h> #include <CGAL/Nef_3/SNC_indexed_items.h> #undef CGAL_NEF_DEBUG #define CGAL_NEF_DEBUG 29 #include <CGAL/Nef_2/debug.h> CGAL_BEGIN_NAMESPACE template <class SNC_structure> class Shell_to_nef_3 { typedef typename SNC_structure::SM_decorator SM_decorator; typedef typename SNC_structure::Vertex_handle Vertex_handle; typedef typename SNC_structure::SVertex_handle SVertex_handle; typedef typename SNC_structure::SHalfedge_handle SHalfedge_handle; typedef typename SNC_structure::SFace_handle SFace_handle; typedef typename SNC_structure::Point_3 Point_3; typedef typename SNC_structure::Plane_3 Plane_3; typedef typename SNC_structure::Sphere_point Sphere_point; typedef typename SNC_structure::Sphere_segment Sphere_segment; typedef typename SNC_structure::Sphere_circle Sphere_circle; typedef typename SNC_structure::Vertex_const_handle Vertex_const_handle; typedef typename SNC_structure::Halfedge_const_handle Halfedge_const_handle; typedef typename SNC_structure::Halffacet_const_handle Halffacet_const_handle; typedef typename SNC_structure::SHalfedge_const_handle SHalfedge_const_handle; typedef typename SNC_structure::SHalfloop_const_handle SHalfloop_const_handle; typedef typename SNC_structure::SFace_const_handle SFace_const_handle; typedef typename SNC_structure::SFace_cycle_const_iterator SFace_cycle_const_iterator; typedef typename SNC_structure::SHalfedge_around_sface_const_circulator SHalfedge_around_sface_const_circulator; /* template<typename Items, typename Polyhedron, typename SNC_structure> class Index_adder { typedef typename SNC_structure::SHalfedge_handle SHalfedge_handle; typedef typename Polyhedron::Halfedge_around_vertex_const_circulator Halfedge_around_vertex_const_circulator; public: Index_adder(Polyhedron& ) {} void set_hash(Halfedge_around_vertex_const_circulator, SHalfedge_handle) {} void resolve_indexes() {} }; template<typename Polyhedron, typename SNC_structure> class Index_adder<CGAL::SNC_indexed_items, Polyhedron, SNC_structure> { typedef typename SNC_structure::SHalfedge_handle SHalfedge_handle; typedef typename Polyhedron::Halfedge_const_handle Halfedge_const_handle; typedef typename Polyhedron::Facet_const_iterator Facet_const_iterator; typedef typename Polyhedron::Halfedge_around_vertex_const_circulator Halfedge_around_vertex_const_circulator; typedef typename Polyhedron::Halfedge_around_facet_const_circulator Halfedge_around_facet_const_circulator; typedef typename CGAL::Unique_hash_map<Halfedge_const_handle, SHalfedge_handle> Hash; Polyhedron& P; Hash hash; public: Index_adder(Polyhedron& P_) : P(P_) {} void set_hash(Halfedge_around_vertex_const_circulator evc, SHalfedge_handle se) { hash[evc] = se; } void resolve_indexes() { Facet_const_iterator fi; for(fi = P.facets_begin(); fi != P.facets_end(); ++fi) { Halfedge_around_facet_const_circulator fc(fi->facet_begin()), end(fc); hash[fc]->set_index(); hash[fc]->twin()->set_index(); hash[fc]->twin()->source()->set_index(); int se = hash[fc]->get_index(); int set = hash[fc]->twin()->get_index(); int sv = hash[fc]->twin()->source()->get_index(); ++fc; CGAL_For_all(fc, end) { hash[fc]->set_index(se); hash[fc]->twin()->set_index(set); hash[fc]->source()->set_index(sv); hash[fc]->twin()->source()->set_index(); sv = hash[fc]->twin()->source()->get_index(); } hash[fc]->source()->set_index(sv); } } }; */ private: SNC_structure& S; public: Shell_to_nef_3(SNC_structure& S_) : S(S_) /* : index_adder(sf)*/ {} void visit(Vertex_const_handle ) {} void visit(Halfedge_const_handle ) {} void visit(Halffacet_const_handle ) {} void visit(SHalfedge_const_handle ) {} void visit(SHalfloop_const_handle ) {} void visit(SFace_const_handle sf) { SFace_cycle_const_iterator sfci = sf->sface_cycles_begin(); if(sfci.is_shalfloop()) return; CGAL_assertion(sfci.is_shalfedge()); Vertex_const_handle pv = sf->center_vertex(); Vertex_handle nv = S.new_vertex(); nv->point() = pv->point(); nv->mark() = true; SM_decorator SM(&*nv); SHalfedge_around_sface_const_circulator pe(sfci), pe_prev(pe), pend(pe); SVertex_handle sv_0 = SM.new_svertex(pe->source()->point()); sv_0->mark() = true; ++pe; SVertex_handle sv_prev = sv_0; do { SVertex_handle sv = SM.new_svertex(pe->source()->point()); sv->mark() = true; SHalfedge_handle e = SM.new_shalfedge_pair(sv_prev, sv); e->circle() = pe_prev->circle(); e->twin()->circle() = pe_prev->twin()->circle(); e->mark() = e->twin()->mark() = true; // index_adder.set_hash(pe_prev, e); sv_prev = sv; pe_prev = pe; ++pe; } while( pe != pend ); SHalfedge_handle e; e = SM.new_shalfedge_pair(sv_prev, sv_0); e->circle() = pe_prev->circle(); e->twin()->circle() = pe_prev->twin()->circle(); e->mark() = e->twin()->mark() = true; // index_adder.set_hash(pe_prev, e); // create faces SFace_handle fext = SM.new_sface(); SM.link_as_face_cycle(e->twin(), fext); fext->mark() = false; SFace_handle fint = SM.new_sface(); SM.link_as_face_cycle(e, fint); fint->mark() = false; SM.check_integrity_and_topological_planarity(); } }; template <class Nef_polyhedron> void shell_to_nef_3(const Nef_polyhedron& N, typename Nef_polyhedron::SFace_const_handle sf, typename Nef_polyhedron::SNC_structure& S) { Shell_to_nef_3<typename Nef_polyhedron::SNC_structure> s2n(S); N.visit_shell_objects(sf, s2n); // s2n.resolve_indexes(); } CGAL_END_NAMESPACE #endif //CGAL_NEF_SHELL_TO_NEF_3_H
#include <CGAL/Gmpz.h> #include <CGAL/Homogeneous.h> #include <CGAL/Nef_polyhedron_3.h> #include <CGAL/Polyhedron_3.h> #include <CGAL/IO/Nef_polyhedron_iostream_3.h> int main() { typedef CGAL::Homogeneous<CGAL::Gmpz> Kernel; typedef CGAL::Nef_polyhedron_3<Kernel> Nef_3; typedef CGAL::Polyhedron_3<Kernel> Polyhedron_3; Nef_3 N; std::cin >> N; Nef_3::Volume_const_iterator c = ++(N.volumes_begin()); Nef_3::SFace_const_handle sf = c->shells_begin(); CGAL::Timer t; t.start(); Polyhedron_3 P; N.convert_inner_shell_to_polyhedron(sf, P); Nef_3 R1(P); t.stop(); std::cerr << "old impl " << t.time() << std::endl; t.reset(); t.start(); Nef_3 R2(N, sf); t.stop(); std::cerr << "new impl " << t.time() << std::endl; CGAL_assertion(R1 == R2); CGAL_assertion(R2.is_valid()); return 0; }
- [cgal-discuss] Coordinate Transformation, Max, 01/17/2008
- Re: [cgal-discuss] Coordinate Transformation, Max, 01/17/2008
- [cgal-discuss][Nef_polehedron]convert_to_Polyhedron ( Polyhedron& P), Max, 01/17/2008
- Re: [cgal-discuss][Nef_polehedron]convert_to_Polyhedron ( Polyhedron& P), Peter Hachenberger, 01/17/2008
- Re: Re: [cgal-discuss][Nef_polyhedron]convert_to_Polyhedron (Polyhedron& P), Max, 01/18/2008
- Re: Re: Re: [cgal-discuss][Nef_polyhedron]convert_to_Polyhedron (Polyhedron& P), Max, 01/18/2008
- Re: Re: Re: [cgal-discuss][Nef_polyhedron]convert_to_Polyhedron (Polyhedron& P), Peter Hachenberger, 01/18/2008
- Re: Re: Re: [cgal-discuss][Nef_polyhedron]convert_to_Polyhedron (Polyhedron& P), Max, 01/18/2008
- [cgal-discuss] Order of including CGAL header files, Max, 01/18/2008
- Re: Re: [cgal-discuss][Nef_polyhedron]convert_to_Polyhedron (Polyhedron& P), Max, 01/18/2008
- Re: [cgal-discuss][Nef_polehedron]convert_to_Polyhedron ( Polyhedron& P), Peter Hachenberger, 01/17/2008
Archive powered by MHonArc 2.6.16.