BWAPI
|
00001 // Copyright (c) 1997-2000 Max-Planck-Institute Saarbruecken (Germany). 00002 // All rights reserved. 00003 // 00004 // This file is part of CGAL (www.cgal.org); you may redistribute it under 00005 // the terms of the Q Public License version 1.0. 00006 // See the file LICENSE.QPL distributed with CGAL. 00007 // 00008 // Licensees holding a valid commercial license may use this file in 00009 // accordance with the commercial license agreement provided with the software. 00010 // 00011 // This file is provided AS IS with NO WARRANTY OF ANY KIND, INCLUDING THE 00012 // WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE. 00013 // 00014 // $URL: svn+ssh://scm.gforge.inria.fr/svn/cgal/branches/CGAL-3.5-branch/Nef_S2/include/CGAL/Nef_S2/SM_triangulator.h $ 00015 // $Id: SM_triangulator.h 46617 2008-11-02 20:34:48Z afabri $ 00016 // 00017 // 00018 // Author(s) : Michael Seel <seel@mpi-sb.mpg.de> 00019 00020 #ifndef CGAL_SM_TRIANGULATOR_H 00021 #define CGAL_SM_TRIANGULATOR_H 00022 00023 #include <CGAL/basic.h> 00024 #include <CGAL/Unique_hash_map.h> 00025 #include <CGAL/Nef_2/Segment_overlay_traits.h> 00026 #include <CGAL/Nef_2/geninfo.h> 00027 #include <CGAL/Nef_S2/SM_decorator.h> 00028 #include <CGAL/Nef_S2/SM_const_decorator.h> 00029 #include <CGAL/Nef_S2/SM_point_locator.h> 00030 #include <CGAL/Nef_S2/SM_io_parser.h> 00031 #include <CGAL/Nef_S2/SM_constrained_triang_traits.h> 00032 00033 #undef CGAL_NEF_DEBUG 00034 #define CGAL_NEF_DEBUG 137 00035 #include <CGAL/Nef_2/debug.h> 00036 00037 #define CGAL_USING(t) typedef typename Base::t t 00038 #ifndef CGAL_USE_LEDA 00039 #define LEDA_MEMORY(t) 00040 #endif 00041 CGAL_BEGIN_NAMESPACE 00042 00043 template <typename Decorator_, typename IT, typename INFO> 00044 struct SM_subdivision { 00045 typedef Decorator_ Triangulator; 00046 typedef typename Decorator_::SVertex_handle Vertex_handle; 00047 typedef typename Decorator_::SHalfedge_handle Halfedge_handle; 00048 typedef typename Decorator_::Sphere_point Point; 00049 typedef typename Decorator_::Sphere_segment Segment; 00050 Triangulator T; 00051 CGAL::Unique_hash_map<IT,INFO>& M; 00052 /* M stores the object that supports the segment that 00053 is input object of the sweep */ 00054 00055 SM_subdivision(Triangulator Ti, 00056 CGAL::Unique_hash_map<IT,INFO>& Mi) : T(Ti), M(Mi) {} 00057 00058 Vertex_handle new_vertex(const Point& p) 00059 { Vertex_handle v = T.new_svertex(p); T.assoc_info(v); 00060 return v; 00061 } 00062 00063 void link_as_target_and_append(Vertex_handle v, Halfedge_handle e) 00064 { T.link_as_target_and_append(v,e); } 00065 00066 Halfedge_handle new_halfedge_pair_at_source(Vertex_handle v) 00067 { Halfedge_handle e = 00068 T.new_shalfedge_pair_at_source(v,Decorator_::BEFORE); 00069 T.assoc_info(e); 00070 return e; 00071 } 00072 00073 void halfedge_below(Vertex_handle v, Halfedge_handle e) const 00074 { T.halfedge_below(v) = e; } 00075 00076 /* the following operation associates segment support with 00077 halfedges, we only update if non-NULL; this prevents 00078 artificial sphere subdivision segments that have NULL 00079 support to overwrite non-NULL support */ 00080 00081 void supporting_segment(Halfedge_handle e, IT it) const 00082 { T.is_forward(e) = true; 00083 if ( ! M[it].empty() ) T.support(e) = M[it]; } 00084 00085 /* the following operation associate segment support with 00086 vertices, we only update if non-NULL; this prevents 00087 artificial segments that have NULL support to overwrite 00088 non-NULL support */ 00089 00090 void trivial_segment(Vertex_handle v, IT it) const 00091 { if ( ! M[it].empty() ) T.support(v) = M[it]; } 00092 00093 void starting_segment(Vertex_handle v, IT it) const 00094 { if ( ! M[it].empty() ) T.support(v) = M[it]; } 00095 00096 void ending_segment(Vertex_handle v, IT it) const 00097 { if ( ! M[it].empty() ) T.support(v) = M[it]; } 00098 00099 void passing_segment(Vertex_handle v, IT it) const 00100 { if ( ! M[it].empty() ) T.support(v) = M[it]; } 00101 00102 00103 }; // SM_subdivision 00104 00105 00106 00107 /*{\Manpage {SM_triangulator}{Decorator_}{Overlay in the sphere}{O}}*/ 00108 00109 template <typename Decorator_> 00110 class SM_triangulator : public Decorator_ { 00111 public: 00112 /*{\Mdefinition An instance |\Mvar| of data type |\Mname| is a 00113 decorator object offering sphere map triangulation calculation.}*/ 00114 00115 typedef Decorator_ Base; 00116 typedef typename Decorator_::Sphere_map Sphere_map; 00117 typedef CGAL::SM_const_decorator<Sphere_map> SM_const_decorator; 00118 typedef SM_const_decorator Explorer; 00119 typedef Decorator_ Decorator; 00120 typedef SM_triangulator<Decorator_> Self; 00121 typedef CGAL::SM_point_locator<SM_const_decorator> SM_point_locator; 00122 00123 typedef typename SM_const_decorator::SVertex_const_handle SVertex_const_handle; 00124 typedef typename SM_const_decorator::SHalfedge_const_handle SHalfedge_const_handle; 00125 typedef typename SM_const_decorator::SHalfloop_const_handle SHalfloop_const_handle; 00126 typedef typename SM_const_decorator::SFace_const_handle SFace_const_handle; 00127 typedef typename SM_const_decorator::SVertex_const_iterator SVertex_const_iterator; 00128 typedef typename SM_const_decorator::SHalfedge_const_iterator SHalfedge_const_iterator; 00129 typedef typename SM_const_decorator::SFace_const_iterator SFace_const_iterator; 00130 00131 typedef typename Base::SVertex_handle SVertex_handle; 00132 typedef typename Base::SHalfedge_handle SHalfedge_handle; 00133 typedef typename Base::SHalfloop_handle SHalfloop_handle; 00134 typedef typename Base::SFace_handle SFace_handle; 00135 typedef typename Base::SVertex_iterator SVertex_iterator; 00136 typedef typename Base::SHalfedge_iterator SHalfedge_iterator; 00137 typedef typename Base::SFace_iterator SFace_iterator; 00138 typedef typename Base::Object_handle Object_handle; 00139 00140 typedef typename Base::SHalfedge_around_svertex_circulator 00141 SHalfedge_around_svertex_circulator; 00142 typedef typename Base::SHalfedge_around_sface_circulator 00143 SHalfedge_around_sface_circulator; 00144 00145 typedef std::pair<SHalfedge_handle,SHalfedge_handle> SHalfedge_pair; 00146 00147 /*{\Mtypes 3}*/ 00148 00149 typedef typename Base::Sphere_kernel Sphere_kernel; 00150 00151 typedef typename Sphere_kernel::Sphere_point Sphere_point; 00152 /*{\Mtypemember the point type of the sphere geometry.}*/ 00153 typedef typename Sphere_kernel::Sphere_segment Sphere_segment; 00154 /*{\Mtypemember the segment type of the sphere geometry.}*/ 00155 typedef typename Sphere_kernel::Sphere_circle Sphere_circle; 00156 /*{\Mtypemember the circle type of the sphere geometry.}*/ 00157 typedef typename Sphere_kernel::Sphere_triangle Sphere_triangle; 00158 /*{\Mtypemember the triangle type of the sphere geometry.}*/ 00159 00160 typedef typename Decorator::Mark Mark; 00161 /*{\Mtypemember the mark of sphere map objects.}*/ 00162 00163 /*{\Mgeneralization Decorator_}*/ 00164 00165 protected: 00166 const Explorer* E_; 00167 const Sphere_kernel& K; 00168 00169 public: 00170 00171 typedef std::list<Sphere_segment> Seg_list; 00172 typedef typename Seg_list::iterator Seg_iterator; 00173 typedef std::pair<Seg_iterator,Seg_iterator> Seg_it_pair; 00174 typedef std::pair<Sphere_segment,Sphere_segment> Seg_pair; 00175 typedef CGAL::Unique_hash_map<Seg_iterator,Object_handle> Seg_map; 00176 00177 // vertex_info stores the origin of vertices 00178 struct vertex_info { 00179 Object_handle o_; 00180 SHalfedge_handle e_; 00181 vertex_info() : o_(),e_() {} 00182 LEDA_MEMORY(vertex_info) 00183 }; 00184 00185 void assoc_info(SVertex_handle v) const 00186 { geninfo<vertex_info>::create(info(v)); } 00187 00188 void discard_info(SVertex_handle v) const 00189 { geninfo<vertex_info>::clear(info(v)); } 00190 00191 vertex_info& ginfo(SVertex_handle v) const 00192 { return geninfo<vertex_info>::access(info(v)); } 00193 00194 Object_handle& support(SVertex_handle v) const 00195 { return ginfo(v).o_; } 00196 00197 SHalfedge_handle& halfedge_below(SVertex_handle v) const 00198 { return ginfo(v).e_; } 00199 00200 // edge_info stores the origin of edges 00201 struct edge_info { 00202 Mark m_left_; Object_handle o_; bool forw_; 00203 00204 edge_info() { m_left_=Mark(); o_=Object_handle(); forw_=false; } 00205 LEDA_MEMORY(edge_info) 00206 }; 00207 00208 void assoc_info(SHalfedge_handle e) const 00209 { geninfo<edge_info>::create(info(e)); 00210 geninfo<edge_info>::create(info(e->twin())); } 00211 00212 void discard_info(SHalfedge_handle e) const 00213 { geninfo<edge_info>::clear(info(e)); 00214 geninfo<edge_info>::clear(info(e->twin())); } 00215 00216 edge_info& ginfo(SHalfedge_handle e) const 00217 { return geninfo<edge_info>::access(info(e)); } 00218 00219 Object_handle& support(SHalfedge_handle e) const 00220 // uedge information we store in the smaller one 00221 { if (&*e < &*(e->twin())) return ginfo(e).o_; 00222 else return ginfo(e->twin()).o_; } 00223 00224 Mark& incident_mark(SHalfedge_handle e) const 00225 // biedge information we store in the edge 00226 { return ginfo(e).m_left_; } 00227 00228 const edge_info& ginfo(SHalfedge_const_handle e) const 00229 { return geninfo<edge_info>::const_access(info(e)); } 00230 const Mark& incident_mark(SHalfedge_const_handle e) const 00231 { return ginfo(e).m_left_; } 00232 00233 bool& is_forward(SHalfedge_handle e) const 00234 // biedge information we store in the edge 00235 { return ginfo(e).forw_; } 00236 00237 void assert_equal_marks(SVertex_handle v1, SVertex_handle v2) const 00238 { CGAL_assertion(v1->mark()==v2->mark()); } 00239 00240 void assert_equal_marks(SHalfedge_handle e1, SHalfedge_handle e2) const 00241 { CGAL_assertion(e1->mark()==e2->mark()); } 00242 00243 Sphere_segment segment(const Explorer* , 00244 SHalfedge_const_handle e) const 00245 { return Sphere_segment( 00246 e->source()->point(),e->twin()->source()->point(),e->circle()); } 00247 00248 Sphere_segment trivial_segment(const Explorer* , 00249 SVertex_const_handle v) const 00250 { Sphere_point p = v->point(); 00251 return Sphere_segment(p,p); } 00252 00253 Seg_pair two_segments(const Explorer* , 00254 SHalfedge_const_handle e) const 00255 // we know that source(e)==target(e) 00256 { return e->circle().split_at(e->source()->point()); } 00257 00258 Seg_pair two_segments(const Explorer* , 00259 SHalfloop_const_handle l) const 00260 { return l->circle().split_at_xy_plane(); } 00261 00262 /*{\Mcreation 6}*/ 00263 SM_triangulator(Sphere_map* M, const Explorer* E, 00264 const Sphere_kernel& G = Sphere_kernel()) : Base(M), E_(E), K(G) {} 00265 /*{\Mcreate |\Mvar| is a triangulator object for the map |M|, 00266 stores the triangulation in |MT|.}*/ 00267 00268 /*{\Moperations 1.1 1}*/ 00269 00270 void triangulate(); 00271 /*{\Mop produces a triangulated sphere map.}*/ 00272 00273 void triangulate_per_hemisphere(SVertex_iterator start, SVertex_iterator end); 00274 00275 template <typename Iterator, typename T> 00276 void partition_to_halfsphere(Iterator start, Iterator end, 00277 Seg_list& L, CGAL::Unique_hash_map<Iterator,T>& M, int pos) const; 00278 00279 void merge_halfsphere_maps(SVertex_handle v1, SVertex_handle v2); 00280 void merge_nodes(SHalfedge_handle e1, SHalfedge_handle e2); 00281 void complete_support(SVertex_iterator v_start, SVertex_iterator v_end, 00282 Mark mohs) const; 00283 00284 void correct_triangle_at(SVertex_handle v) 00285 { CGAL_NEF_TRACEN("correct_triangle_at "<<PH(v)); 00286 if ( !has_outdeg_two(v) ) return; 00287 SHalfedge_handle e = first_out_edge(v); 00288 CGAL_assertion(e->snext()->snext()->snext()==e); 00289 flip_diagonal(e->snext()); 00290 } 00291 00292 void dump(std::ostream& os = std::cerr) const 00293 { SM_io_parser<Explorer>::dump(E_,os); 00294 SM_io_parser<Base>::dump(*this,os); } 00295 00296 Sphere_triangle incident_triangle(SHalfedge_handle e) const 00297 { SHalfedge_handle en(e->snext()), enn(en->snext()); 00298 CGAL_assertion(enn->snext()==e); 00299 return Sphere_triangle(e->source()->point(), 00300 en->source()->point(), 00301 enn->source()->point(), 00302 e->circle(), 00303 en->circle(), 00304 enn->circle()); 00305 } 00306 00307 Sphere_triangle incident_triangle(SHalfedge_const_handle e) const 00308 { SHalfedge_const_handle en(e->snext()), enn(en->snext()); 00309 CGAL_assertion(enn->snext()==e); 00310 return Sphere_triangle(e->source()->point(), 00311 en->source()->point(), 00312 enn->source()->point(), 00313 e->circle(), 00314 en->circle(), 00315 enn->circle()); 00316 } 00317 00318 void discard_info() 00319 { 00320 SVertex_iterator v; 00321 SHalfedge_iterator e; 00322 CGAL_forall_svertices(v,*this) discard_info(v); 00323 CGAL_forall_shalfedges(e,*this) discard_info(e); 00324 } 00325 00326 00327 }; // SM_triangulator<Decorator_> 00328 00329 00330 template <typename Decorator_> 00331 void SM_triangulator<Decorator_>::triangulate() 00332 { CGAL_NEF_TRACEN("triangulate"); 00333 // first create sphere segments from isoverts, edges, loops 00334 Seg_list L; 00335 Seg_map From; 00336 SVertex_const_iterator v; 00337 CGAL_forall_svertices(v,*E_) { 00338 if ( !E_->is_isolated(v) ) continue; 00339 L.push_back(trivial_segment(E_,v)); 00340 From[--L.end()] = make_object(v); 00341 } 00342 SHalfedge_const_iterator e; 00343 CGAL_forall_sedges(e,*E_) { 00344 if ( e->source() == e->twin()->source() ) { 00345 Seg_pair p = two_segments(E_,e); 00346 L.push_back(p.first); L.push_back(p.second); 00347 From[--L.end()] = From[--(--L.end())] = make_object(e); 00348 } else { 00349 L.push_back(segment(E_,e)); 00350 From[--L.end()] = make_object(e); 00351 } 00352 } 00353 if ( E_->has_shalfloop() ) { 00354 Seg_pair p = two_segments(E_,E_->shalfloop()); 00355 L.push_back(p.first); L.push_back(p.second); 00356 From[--L.end()] = From[--(--L.end())] = make_object(E_->shalfloop()); 00357 } 00358 00359 // partition segments from L to positive and negative hemisphere 00360 Seg_list L_pos,L_neg; 00361 partition_to_halfsphere(L.begin(), L.end(), L_pos, From, +1); 00362 partition_to_halfsphere(L.begin(), L.end(), L_neg, From, -1); 00363 00364 // typename Seg_list::iterator it; 00365 // std::cerr << "L_pos" << std::endl; 00366 // CGAL_forall_iterators(it,L_pos) std::cerr << *it << std::endl; 00367 // std::cerr << "L_neg" << std::endl; 00368 // CGAL_forall_iterators(it,L_neg) std::cerr << *it << std::endl; 00369 00370 // sweep the hemispheres to create two half sphere maps 00371 typedef SM_subdivision<Self,Seg_iterator,Object_handle> SM_output; 00372 typedef typename Sphere_kernel::Positive_halfsphere_geometry PH_geometry; 00373 typedef CGAL::Segment_overlay_traits< 00374 Seg_iterator, SM_output, PH_geometry> PHS_traits; 00375 typedef CGAL::generic_sweep<PHS_traits> Positive_halfsphere_sweep; 00376 00377 typedef typename Sphere_kernel::Negative_halfsphere_geometry NH_geometry; 00378 typedef CGAL::Segment_overlay_traits< 00379 Seg_iterator, SM_output, NH_geometry> NHS_traits; 00380 typedef CGAL::generic_sweep<NHS_traits> Negative_halfsphere_sweep; 00381 00382 SVertex_handle v_sep; 00383 SHalfedge_handle e_sep; 00384 SM_output O(*this,From); 00385 00386 typedef typename PHS_traits::INPUT Input_range; 00387 Positive_halfsphere_sweep SP( 00388 Input_range(L_pos.begin(),L_pos.end()),O, 00389 PH_geometry()); 00390 SP.sweep(); 00391 v_sep=--this->svertices_end(); e_sep=--this->shalfedges_end(); 00392 00393 Negative_halfsphere_sweep SM( 00394 Input_range(L_neg.begin(),L_neg.end()),O, 00395 NH_geometry()); 00396 SM.sweep(); 00397 ++v_sep; ++e_sep; 00398 // now two CCs of sphere graph are calculated 00399 // v_sep = first vertex of CC in negative x-sphere 00400 // e_sep = first edge of CC in negative x-sphere 00401 00402 SHalfedge_iterator u; 00403 CGAL_forall_sedges(u,*this) { 00404 Sphere_segment s(u->source()->point(),u->twin()->source()->point()); 00405 u->circle() = s.sphere_circle(); 00406 u->twin()->circle() = s.sphere_circle().opposite(); 00407 } 00408 00409 Mark lower, upper; 00410 SM_point_locator PL(E_->sphere_map()); 00411 PL.marks_of_halfspheres(lower,upper); 00412 complete_support(this->svertices_begin(), v_sep, lower); 00413 complete_support(v_sep, this->svertices_end(), upper); 00414 00415 /* 00416 CGAL_forall_sedges(u,*this) { 00417 std::cerr << point(source(u)) << "->" << point(target(u)) << std::endl; 00418 } 00419 */ 00420 00421 // triangulate per hemisphere 00422 typedef SM_constrained_triang_traits<Self,PH_geometry> PCT_traits; 00423 typedef CGAL::generic_sweep<PCT_traits> Positive_halfsphere_ct_sweep; 00424 typedef SM_constrained_triang_traits<Self,NH_geometry> NCT_traits; 00425 typedef CGAL::generic_sweep<NCT_traits> Negative_halfsphere_ct_sweep; 00426 typedef std::pair<SVertex_iterator,SVertex_iterator> SVertex_pair; 00427 00428 SVertex_pair vpp(this->svertices_begin(),v_sep); 00429 Positive_halfsphere_ct_sweep PCTS(vpp, *this, 00430 PH_geometry()); 00431 PCTS.sweep(); 00432 SVertex_pair vpn(v_sep,this->svertices_end()); 00433 Negative_halfsphere_ct_sweep NCTS(vpn, *this, 00434 NH_geometry()); 00435 NCTS.sweep(); 00436 00437 /* 00438 std::cerr << std::endl; 00439 CGAL_forall_sedges(u,*this) { 00440 std::cerr << point(source(u)) << "->" << point(target(u)) << std::endl; 00441 } 00442 */ 00443 00444 /* Note the we divide the world along the xy equator and 00445 split the equator at y- and y+. We treat the halfcircle 00446 at x+ as if perturbed slightly up. This makes triangles 00447 that have y- or y+ as a vertex degenerate. if such triangles 00448 appear we repair it by flipping the edge opposite to the 00449 vertex y-(y+). 00450 */ 00451 00452 correct_triangle_at(this->svertices_begin()); 00453 correct_triangle_at(--SVertex_iterator(v_sep)); 00454 correct_triangle_at(v_sep); 00455 correct_triangle_at(--this->svertices_end()); 00456 00457 CGAL_forall_sedges(u,*this) { 00458 Sphere_segment s(u->source()->point(),u->twin()->source()->point()); 00459 u->circle() = s.sphere_circle(); 00460 u->twin()->circle() = s.sphere_circle().opposite(); 00461 } 00462 00463 // merge the hemisphere maps into one sphere map 00464 merge_halfsphere_maps(this->svertices_begin(),v_sep); 00465 this->check_integrity_and_topological_planarity(false); 00466 } 00467 00468 00469 template <typename Decorator_> 00470 template <typename Iterator, typename T> 00471 void SM_triangulator<Decorator_>:: 00472 partition_to_halfsphere(Iterator start, Iterator beyond, Seg_list& L, 00473 CGAL::Unique_hash_map<Iterator,T>& M, int pos) const 00474 { CGAL_NEF_TRACEN("partition_to_halfsphere "); 00475 CGAL_assertion(pos!=0); 00476 bool add_cross = true; 00477 Sphere_segment s1,s2; 00478 Sphere_circle xycircle(0,0,pos); 00479 while ( start != beyond ) { 00480 if(start->source().hz() * pos > 0 || start->target().hz() * pos > 0) 00481 add_cross = false; 00482 int i = start->intersection(xycircle,s1,s2); 00483 if (i>1) { L.push_back(s2); M[--L.end()] = M[start]; } 00484 if (i>0) { L.push_back(s1); M[--L.end()] = M[start]; } 00485 ++start; 00486 } 00487 // now all segments are split into halfspheres 00488 // we still have to: 00489 // - split segments containing our special poles y^-, y^+ 00490 // - split halfcircles 00491 // - add four equator segments 00492 Sphere_point S(0,-1,0),N(0,1,0); 00493 Sphere_circle yzcircle(1,0,0); 00494 typename Seg_list::iterator it, itl; 00495 00496 bool part_in_hemisphere(false); 00497 CGAL_forall_iterators(it,L) { CGAL_NEF_TRACEN(" "<<*it); 00498 if ( equal_as_sets(it->sphere_circle(),xycircle) ) { 00499 CGAL_NEF_TRACEN(" splitting xy seg "<<*it); 00500 int n1 = it->intersection(yzcircle,s1,s2); 00501 if (n1 > 1 && !s2.is_degenerate()) 00502 { M[ L.insert(it,s2) ] = M[it]; } 00503 if (n1 > 0 && !s1.is_degenerate()) 00504 { M[ L.insert(it,s1) ] = M[it]; } 00505 int n2 = it->intersection(yzcircle.opposite(),s1,s2); 00506 if (n2 > 1 && !s2.is_degenerate()) 00507 { M[ L.insert(it,s2) ] = M[it]; } 00508 if (n2 > 0 && !s1.is_degenerate()) 00509 { M[ L.insert(it,s1) ] = M[it]; } 00510 itl = it; --it; L.erase(itl); M[itl] = T(); 00511 // at least one item was appended 00512 } else { 00513 part_in_hemisphere = true; 00514 } 00515 } 00516 CGAL_forall_iterators(it,L) { 00517 if ( it->is_halfcircle() ) { 00518 CGAL_NEF_TRACEN(" splitting halfcircle "<<*it); 00519 Sphere_segment s1,s2; 00520 it->split_halfcircle(s1,s2); 00521 *it = s2; 00522 M[ L.insert(it,s1) ] = M[it]; 00523 } 00524 } 00525 // append 4 xy-equator segments: 00526 Sphere_segment sp(S,N,xycircle); 00527 Sphere_segment sm(S,N,xycircle.opposite()); 00528 Sphere_segment s[4]; 00529 sp.split_halfcircle(s[0],s[1]); 00530 sm.split_halfcircle(s[2],s[3]); 00531 L.insert(L.end(),s,s+4); 00532 /* if no segment is covering the interior of the hemisphere 00533 we have to add a trivial segment to allow for a correct 00534 triangulation */ 00535 if ( !part_in_hemisphere || add_cross) { 00536 Sphere_point p(0,0,pos); 00537 Sphere_circle c(1,0,0); 00538 L.push_back(Sphere_segment(p,p,c)); 00539 } 00540 } 00541 00542 template <typename Decorator_> 00543 void SM_triangulator<Decorator_>:: 00544 merge_nodes(SHalfedge_handle e1, SHalfedge_handle e2) 00545 { 00546 SVertex_handle v1 = e1->source(), v2 = e2->twin()->source(); 00547 CGAL_NEF_TRACEN("merge_nodes "<<PH(v1)<<PH(v2)); 00548 CGAL_assertion(v1->point()==v2->point()); 00549 SHalfedge_handle ep1 = e1->sprev(), en2 = e2->snext(); 00550 SHalfedge_around_svertex_circulator eav(out_edges(v2)),ee(eav); 00551 CGAL_For_all(eav,ee) { set_source(eav,v1); } 00552 link_as_prev_next_pair(e2,e1); 00553 link_as_prev_next_pair(ep1,en2); 00554 assert_equal_marks(v1,v2); 00555 discard_info(v2); 00556 delete_vertex_only(v2); 00557 } 00558 00559 00560 template <typename Decorator_> 00561 void SM_triangulator<Decorator_>:: 00562 merge_halfsphere_maps(SVertex_handle v1, SVertex_handle v2) 00563 { CGAL_NEF_TRACEN("merging halfspheres "<<PH(v1)<<PH(v2)); 00564 CGAL_assertion(v1->point()==v2->point()); 00565 std::list<SHalfedge_pair> L_equator; 00566 SHalfedge_around_sface_circulator 00567 ep(last_out_edge(v1)), en(first_out_edge(v2)->twin()); 00568 do { 00569 L_equator.push_back(SHalfedge_pair(ep,en)); 00570 merge_nodes(ep,en); ++ep; --en; 00571 } while ( ep->source() != v1 ); 00572 00573 typename std::list<SHalfedge_pair>::iterator it; 00574 CGAL_forall_iterators(it,L_equator) { 00575 SHalfedge_handle e1 = it->first, e2 = it->second; 00576 SHalfedge_handle e1t = e1->twin(), e2t = e2->twin(); 00577 CGAL_NEF_TRACEV(PH(e1));CGAL_NEF_TRACEV(PH(e2)); 00578 SHalfedge_handle e2tp = e2t->sprev(); 00579 SHalfedge_handle e2tn = e2t->snext(); 00580 link_as_prev_next_pair(e2tp,e1); 00581 link_as_prev_next_pair(e1,e2tn); 00582 SFace_handle f = e2t->incident_sface(); 00583 if ( is_sm_boundary_object(e2t) ) 00584 { undo_sm_boundary_object(e2t,f); store_sm_boundary_object(e1,f); } 00585 set_face(e1,f); 00586 if ( e2 == first_out_edge(e2->source()) ) 00587 set_first_out_edge(e2->source(),e1t); 00588 discard_info(e2); 00589 delete_edge_pair_only(e2); 00590 } 00591 } 00592 00593 template <typename Decorator_> 00594 void SM_triangulator<Decorator_>:: 00595 complete_support(SVertex_iterator v_start, SVertex_iterator v_end, 00596 Mark mohs) const 00597 { CGAL_NEF_TRACEN("complete_support"); 00598 Mark m_buffer(mohs); 00599 for (SVertex_iterator v = v_start; v != v_end; ++v) { 00600 CGAL_NEF_TRACEN(" vertex = "<<PH(v)); 00601 SHalfedge_handle e_below = halfedge_below(v); 00602 if ( v != v_start ) { 00603 if ( e_below != SHalfedge_handle() ) { 00604 m_buffer = incident_mark(e_below); 00605 } else { // e_below does not exist 00606 /* this is only the case for a vertex v on the final equatorial 00607 halfcircle; there we take the mark from an inedge edge into v */ 00608 // CGAL_assertion( point(v).z() == 0 && 00609 // ( pos > 0 ? (point(v).x() >= 0) : (point(v).x()<=0)) ); 00610 m_buffer = incident_mark(first_out_edge(v)->sprev()); 00611 } 00612 } 00613 CGAL_NEF_TRACEN(" face mark below "<<m_buffer); 00614 00615 Object_handle o = support(v); 00616 SVertex_const_handle vs; 00617 SHalfedge_const_handle es; 00618 SHalfloop_const_handle ls; 00619 if ( o.empty() ) { v->mark() = m_buffer; } 00620 else if ( CGAL::assign(vs,o) ) { v->mark() = vs->mark(); } 00621 else if ( CGAL::assign(es,o) ) { 00622 if ( es->source()->point() == v->point() ) 00623 { v->mark() = es->source()->mark(); } 00624 else if ( es->twin()->source()->point() == v->point() ) 00625 { v->mark() = es->twin()->source()->mark(); } 00626 else { v->mark() = es->mark(); } 00627 } 00628 else if ( CGAL::assign(ls,o) ) { v->mark() = ls->mark(); } 00629 else CGAL_error_msg("damn wrong support."); 00630 CGAL_NEF_TRACEN(" face mark at "<<v->mark()); 00631 00632 if ( is_isolated(v) ) continue; 00633 00634 SHalfedge_around_svertex_circulator e(first_out_edge(v)), hend(e); 00635 CGAL_For_all(e,hend) { 00636 CGAL_NEF_TRACEN(" edge "<<PH(e)); 00637 if ( !is_forward(e) ) break; 00638 if ( ! support(e).empty() ) { 00639 SHalfedge_const_handle ei; 00640 if ( CGAL::assign(ei,support(e)) ) { 00641 if ( ei->circle() != e->circle() ) { ei = ei->twin(); } 00642 CGAL_assertion( ei->circle() == e->circle() ); 00643 CGAL_NEF_TRACEN(" supporting edge "<<PH(ei)); 00644 incident_mark(e->twin()) = ei->twin()->incident_sface()->mark(); 00645 e->mark() = ei->mark(); 00646 incident_mark(e) = m_buffer = ei->incident_sface()->mark(); 00647 } 00648 SHalfloop_const_handle li; 00649 if ( CGAL::assign(li,support(e)) ) { 00650 if ( li->circle() != e->circle() ) { li = li->twin(); } 00651 CGAL_assertion( li->circle() == e->circle() ); 00652 CGAL_NEF_TRACEN(" supporting loop "<<PH(li)); 00653 incident_mark(e->twin()) = li->twin()->incident_sface()->mark(); 00654 e->mark() = li->mark(); 00655 incident_mark(e) = m_buffer = li->incident_sface()->mark(); 00656 } 00657 } else { CGAL_NEF_TRACEN(" support from face below "); 00658 incident_mark(e->twin()) = e->mark() = 00659 incident_mark(e) = m_buffer; 00660 } 00661 CGAL_NEF_TRACEN(" new face mark "<<m_buffer); 00662 00663 } // CGAL_For_all(e,hend) 00664 00665 CGAL_NEF_TRACEN(" mark of "<<PH(v)); 00666 } 00667 00668 } 00669 00670 00671 00672 00673 CGAL_END_NAMESPACE 00674 #undef CGAL_USING 00675 #endif //CGAL_SM_TRIANGULATOR_H 00676 00677