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_2/include/CGAL/Nef_2/Constrained_triang_traits.h $ 00015 // $Id: Constrained_triang_traits.h 40822 2007-11-07 16:51:18Z ameyer $ 00016 // 00017 // 00018 // Author(s) : Michael Seel <seel@mpi-sb.mpg.de> 00019 #ifndef CGAL_PM_CONSTR_TRIANG_TRAITS_H 00020 #define CGAL_PM_CONSTR_TRIANG_TRAITS_H 00021 00022 #include <CGAL/basic.h> 00023 #include <CGAL/Unique_hash_map.h> 00024 #include <CGAL/generic_sweep.h> 00025 #include <CGAL/Nef_2/PM_checker.h> 00026 #include <cstdlib> 00027 #include <string> 00028 #include <map> 00029 #include <set> 00030 #undef CGAL_NEF_DEBUG 00031 #define CGAL_NEF_DEBUG 19 00032 #include <CGAL/Nef_2/debug.h> 00033 00034 #if defined(BOOST_MSVC) 00035 # pragma warning(push) 00036 # pragma warning(disable:4355) // complaint about using 'this' to 00037 #endif // initialize a member 00038 00039 CGAL_BEGIN_NAMESPACE 00040 00041 struct Do_nothing { 00042 Do_nothing() {} 00043 template <typename ARG> 00044 void operator()(ARG&) const {} 00045 }; 00046 00047 00048 template <typename PMDEC, typename GEOM, 00049 typename NEWEDGE = Do_nothing> 00050 class Constrained_triang_traits : public PMDEC { 00051 public: 00052 typedef Constrained_triang_traits<PMDEC,GEOM,NEWEDGE> Self; 00053 typedef PMDEC Base; 00054 00055 // the types interfacing the sweep: 00056 typedef NEWEDGE INPUT; 00057 typedef typename PMDEC::Plane_map OUTPUT; 00058 typedef GEOM GEOMETRY; 00059 00060 typedef typename GEOM::Point_2 Point; 00061 typedef typename GEOM::Segment_2 Segment; 00062 typedef typename GEOM::Direction_2 Direction; 00063 00064 typedef typename Base::Halfedge_handle Halfedge_handle; 00065 typedef typename Base::Vertex_handle Vertex_handle; 00066 typedef typename Base::Face_handle Face_handle; 00067 typedef typename Base::Halfedge_iterator Halfedge_iterator; 00068 typedef typename Base::Vertex_iterator Vertex_iterator; 00069 typedef typename Base::Face_iterator Face_iterator; 00070 typedef typename Base::Halfedge_base Halfedge_base; 00071 typedef typename Base::Halfedge_around_vertex_circulator 00072 Halfedge_around_vertex_circulator; 00073 00074 class lt_edges_in_sweepline : public PMDEC 00075 { const Point& p; 00076 const Halfedge_handle& e_bottom; 00077 const Halfedge_handle& e_top; 00078 const GEOMETRY& K; 00079 public: 00080 lt_edges_in_sweepline(const Point& pi, 00081 const Halfedge_handle& e1, const Halfedge_handle& e2, 00082 const PMDEC& D, const GEOMETRY& k) : 00083 PMDEC(D), p(pi), e_bottom(e1), e_top(e2), K(k) {} 00084 00085 lt_edges_in_sweepline(const lt_edges_in_sweepline& lt) : 00086 PMDEC(lt), p(lt.p), e_bottom(lt.e_bottom), e_top(lt.e_top), K(lt.K) {} 00087 00088 Segment seg(const Halfedge_handle& e) const 00089 { return K.construct_segment(point(source(e)),point(target(e))); } 00090 00091 int orientation(Halfedge_handle e, const Point& p) const 00092 { return K.orientation(point(source(e)),point(target(e)),p); } 00093 00094 bool operator()(const Halfedge_handle& e1, const Halfedge_handle& e2) const 00095 { // Precondition: 00096 // [[p]] is identical to the source of either [[e1]] or [[e2]]. 00097 if (e1 == e_bottom || e2 == e_top) return true; 00098 if (e2 == e_bottom || e1 == e_top) return false; 00099 if ( e1 == e2 ) return 0; 00100 int s = 0; 00101 if ( p == point(source(e1)) ) s = orientation(e2,p); 00102 else if ( p == point(source(e2)) ) s = - orientation(e1,p); 00103 else CGAL_error_msg("compare error in sweep."); 00104 if ( s || source(e1) == target(e1) || source(e2) == target(e2) ) 00105 return ( s < 0 ); 00106 s = orientation(e2,point(target(e1))); 00107 if (s==0) CGAL_error_msg("parallel edges not allowed."); 00108 return ( s < 0 ); 00109 } 00110 00111 00112 }; // lt_edges_in_sweepline 00113 00114 class lt_pnts_xy : public PMDEC 00115 { const GEOMETRY& K; 00116 public: 00117 lt_pnts_xy(const PMDEC& D, const GEOMETRY& k) : PMDEC(D), K(k) {} 00118 lt_pnts_xy(const lt_pnts_xy& lt) : PMDEC(lt), K(lt.K) {} 00119 int operator()(const Vertex_handle& v1, const Vertex_handle& v2) const 00120 { return K.compare_xy(point(v1),point(v2)) < 0; } 00121 }; // lt_pnts_xy 00122 00123 00124 typedef std::map<Halfedge_handle, Halfedge_handle, lt_edges_in_sweepline> 00125 Sweep_status_structure; 00126 typedef typename Sweep_status_structure::iterator ss_iterator; 00127 typedef typename Sweep_status_structure::value_type ss_pair; 00128 typedef std::set<Vertex_iterator,lt_pnts_xy> Event_Q; 00129 typedef typename Event_Q::const_iterator event_iterator; 00130 00131 const GEOMETRY& K; 00132 Event_Q event_Q; 00133 event_iterator event_it; 00134 Vertex_handle event; 00135 Point p_sweep; 00136 Sweep_status_structure SL; 00137 CGAL::Unique_hash_map<Halfedge_handle,ss_iterator> SLItem; 00138 const NEWEDGE& Treat_new_edge; 00139 Halfedge_handle e_low,e_high; // framing edges ! 00140 Halfedge_handle e_search; 00141 00142 Constrained_triang_traits(const INPUT& in, OUTPUT& out, const GEOMETRY& k) 00143 : Base(out), K(k), event_Q(lt_pnts_xy(*this,K)), 00144 SL(lt_edges_in_sweepline(p_sweep,e_low,e_high,*this,K)), 00145 SLItem(SL.end()), Treat_new_edge(in) 00146 { CGAL_NEF_TRACEN("Constrained Triangulation Sweep"); } 00147 00148 00149 Halfedge_handle new_bi_edge(Vertex_handle v1, Vertex_handle v2) 00150 { // appended at v1 and v2 adj list 00151 Halfedge_handle e = Base::new_halfedge_pair(v1,v2); 00152 Treat_new_edge(e); 00153 return e; 00154 } 00155 00156 Halfedge_handle new_bi_edge(Halfedge_handle e_bf, Halfedge_handle e_af) 00157 { // ccw before e_bf and after e_af 00158 Halfedge_handle e = Base::new_halfedge_pair(e_bf,e_af,Halfedge_base(), 00159 Base::BEFORE, Base::AFTER); 00160 Treat_new_edge(e); 00161 return e; 00162 } 00163 00164 Halfedge_handle new_bi_edge(Vertex_handle v, Halfedge_handle e_bf) 00165 { // appended at v's adj list and before e_bf 00166 Halfedge_handle e = Base::new_halfedge_pair(v,e_bf,Halfedge_base(), 00167 Base::BEFORE); 00168 Treat_new_edge(e); 00169 return e; 00170 } 00171 00172 Segment seg(Halfedge_handle e) const 00173 { return K.construct_segment(point(source(e)),point(target(e))); } 00174 00175 Direction dir(Halfedge_handle e) const 00176 { return K.construct_direction(point(source(e)),point(target(e))); } 00177 00178 bool is_forward(Halfedge_handle e) const 00179 { return K.compare_xy(point(source(e)),point(target(e))) < 0; } 00180 00181 00182 00183 00184 bool edge_is_visible_from(Vertex_handle v, Halfedge_handle e) 00185 { 00186 Point p = point(v); 00187 Point p1 = point(source(e)); 00188 Point p2 = point(target(e)); 00189 return ( K.orientation(p1,p2,p)>0 ); // left_turn 00190 } 00191 00192 void triangulate_up(Halfedge_handle& e_apex) 00193 { 00194 CGAL_NEF_TRACEN("triangulate_up "<<seg(e_apex)); 00195 Vertex_handle v_apex = source(e_apex); 00196 while (true) { 00197 Halfedge_handle e_vis = previous(twin(e_apex)); 00198 bool in_sweep_line = (SLItem[e_vis] != SL.end()); 00199 bool not_visible = !edge_is_visible_from(v_apex,e_vis); 00200 CGAL_NEF_TRACEN(" checking "<<in_sweep_line<<not_visible<<" "<<seg(e_vis)); 00201 if ( in_sweep_line || not_visible) { 00202 CGAL_NEF_TRACEN(" STOP"); return; 00203 } 00204 Halfedge_handle e_back = new_bi_edge(e_apex,e_vis); 00205 if ( !is_forward(e_vis) ) make_first_out_edge(twin(e_back)); 00206 e_apex = e_back; 00207 CGAL_NEF_TRACEN(" produced " << seg(e_apex)); 00208 } 00209 } 00210 00211 void triangulate_down(Halfedge_handle& e_apex) 00212 { 00213 CGAL_NEF_TRACEN("triangulate_down "<<seg(e_apex)); 00214 Vertex_handle v_apex = source(e_apex); 00215 while (true) { 00216 Halfedge_handle e_vis = next(e_apex); 00217 bool in_sweep_line = (SLItem[e_vis] != SL.end()); 00218 bool not_visible = !edge_is_visible_from(v_apex,e_vis); 00219 CGAL_NEF_TRACEN(" checking "<<in_sweep_line<<not_visible<<" "<<seg(e_vis)); 00220 if ( in_sweep_line || not_visible) { 00221 CGAL_NEF_TRACEN(" STOP"); return; 00222 } 00223 Halfedge_handle e_vis_rev = twin(e_vis); 00224 Halfedge_handle e_forw = new_bi_edge(e_vis_rev,e_apex); 00225 e_apex = twin(e_forw); 00226 CGAL_NEF_TRACEN(" produced " << seg(e_apex)); 00227 } 00228 } 00229 00230 void triangulate_between(Halfedge_handle e_upper, Halfedge_handle e_lower) 00231 { 00232 // we triangulate the interior of the whole chain between 00233 // target(e_upper) and target(e_lower) 00234 CGAL_assertion(source(e_upper)==source(e_lower)); 00235 CGAL_NEF_TRACE("triangulate_between\n "<<seg(e_upper)); 00236 CGAL_NEF_TRACEN("\n "<<seg(e_lower)); 00237 Halfedge_handle e_end = twin(e_lower); 00238 while (true) { 00239 Halfedge_handle e_vis = next(e_upper); 00240 Halfedge_handle en_vis = next(e_vis); 00241 CGAL_NEF_TRACEN(" working on base e_vis " << seg(e_vis)); 00242 CGAL_NEF_TRACEN(" next is " << seg(en_vis)); 00243 if (en_vis == e_end) return; 00244 e_upper = twin(new_bi_edge(twin(e_vis),e_upper)); 00245 CGAL_NEF_TRACEN(" produced " << seg(e_upper)); 00246 } 00247 } 00248 00249 void process_event() 00250 { 00251 CGAL_NEF_TRACEN("\nPROCESS_EVENT " << p_sweep); 00252 Halfedge_handle e, ep, eb_low, eb_high, e_end; 00253 if ( !is_isolated(event) ) { 00254 e = last_out_edge(event); 00255 ep = first_out_edge(event); 00256 } 00257 ss_iterator sit_pred, sit; 00258 /* PRECONDITION: 00259 only ingoing => e is lowest in ingoing bundle 00260 only outgoing => e is highest in outgoing bundle 00261 ingoing and outgoing => e is lowest in ingoing bundle */ 00262 eb_high = e_end = ep; 00263 eb_low = e; 00264 CGAL_NEF_TRACEN("determining handle in SL"); 00265 if ( e != Halfedge_handle() ) { 00266 point(target(e_search)) = p_sweep; // degenerate loop edge 00267 sit_pred = SLItem[e]; 00268 if ( sit_pred != SL.end()) sit = --sit_pred; 00269 else sit = sit_pred = --SL.upper_bound(e_search); 00270 } else { // event is isolated vertex 00271 point(target(e_search)) = p_sweep; // degenerate loop edge 00272 sit_pred = --SL.upper_bound(e_search); 00273 } 00274 00275 bool ending_edges(0), starting_edges(0); 00276 while ( e != Halfedge_handle() ) { // walk adjacency list clockwise 00277 if ( SLItem[e] != SL.end() ) 00278 { 00279 CGAL_NEF_TRACEN("ending " << seg(e)); 00280 if (ending_edges) triangulate_between(e,cyclic_adj_succ(e)); 00281 ending_edges = true; 00282 SL.erase(SLItem[e]); 00283 link_bi_edge_to(e,SL.end()); 00284 // not in SL anymore 00285 } 00286 00287 else 00288 { 00289 CGAL_NEF_TRACEN("starting "<<seg(e)); 00290 sit = SL.insert(sit,ss_pair(e,e)); 00291 link_bi_edge_to(e,sit); 00292 if ( !starting_edges ) eb_high = cyclic_adj_succ(e); 00293 starting_edges = true; 00294 } 00295 00296 if (e == e_end) break; 00297 e = cyclic_adj_pred(e); 00298 } 00299 if (!ending_edges) 00300 { 00301 Halfedge_handle e_vis = sit_pred->second; 00302 Halfedge_handle e_vis_n = cyclic_adj_succ(e_vis); 00303 eb_low = eb_high = new_bi_edge(event,e_vis_n); 00304 CGAL_NEF_TRACEN(" producing link "<<seg(eb_low)<<"\n before "<<seg(e_vis_n)); 00305 } 00306 00307 00308 00309 00310 triangulate_up(eb_high); 00311 triangulate_down(eb_low); 00312 sit_pred->second = eb_low; 00313 } 00314 00315 bool event_exists() 00316 { if ( event_it != event_Q.end() ) { 00317 // event is set at end of loop and in init 00318 event = *event_it; 00319 p_sweep = point(event); 00320 return true; 00321 } 00322 return false; 00323 } 00324 00325 void procede_to_next_event() 00326 { ++event_it; } 00327 00328 void link_bi_edge_to(Halfedge_handle e, ss_iterator sit) { 00329 SLItem[e] = SLItem[twin(e)] = sit; 00330 } 00331 00332 void initialize_structures() 00333 { 00334 CGAL_NEF_TRACEN("initialize_structures "); 00335 00336 for ( event=this->vertices_begin(); event != this->vertices_end(); ++event ) 00337 event_Q.insert(event); // sorted order of vertices 00338 00339 event_it = event_Q.begin(); 00340 if ( event_Q.empty() ) return; 00341 event = *event_it; 00342 p_sweep = point(event); 00343 if ( !is_isolated(event) ) { 00344 Halfedge_around_vertex_circulator 00345 e(first_out_edge(event)), eend(e); 00346 CGAL_For_all(e,eend) { 00347 CGAL_NEF_TRACEN("init with "<<PE(e)); 00348 ss_iterator sit = SL.insert(ss_pair(e,e)).first; 00349 link_bi_edge_to(e,sit); 00350 } 00351 } 00352 00353 00354 Vertex_handle v_tmp = this->new_vertex(); point(v_tmp) = Point(); 00355 e_high = Base::new_halfedge_pair(event,v_tmp); 00356 e_low = Base::new_halfedge_pair(event,v_tmp); 00357 // this are two symbolic edges just accessed as sentinels 00358 // they carry no geometric information 00359 e_search = Base::new_halfedge_pair(v_tmp,v_tmp); 00360 // this is just a loop used for searches in SL 00361 00362 ss_iterator sit_high = SL.insert(ss_pair(e_high,e_high)).first; 00363 ss_iterator sit_low = SL.insert(ss_pair(e_low,e_low)).first; 00364 // inserting sentinels into SL 00365 link_bi_edge_to(e_high, sit_high); 00366 link_bi_edge_to(e_low , sit_low); 00367 // we mark them being in the sweepline, which they will never leave 00368 00369 00370 // we move to the second vertex: 00371 procede_to_next_event(); 00372 event_exists(); // sets p_sweep for check invariants 00373 CGAL_NEF_TRACEN("EOF initialization"); 00374 } 00375 00376 void complete_structures() 00377 { 00378 if (e_low != Halfedge_handle()) { 00379 delete_vertex(target(e_search)); 00380 } // removing sentinels and e_search 00381 } 00382 00383 00384 void check_ccw_local_embedding() const 00385 { PM_checker<PMDEC,GEOM> C(*this,K); 00386 C.check_order_preserving_embedding(event); 00387 } 00388 00389 void check_invariants() 00390 { 00391 #ifdef CGAL_CHECK_EXPENSIVE 00392 if ( event_it == event_Q.end() ) return; 00393 check_ccw_local_embedding(); 00394 #endif 00395 } 00396 00397 void check_final() 00398 { 00399 #ifdef CGAL_CHECK_EXPENSIVE 00400 PM_checker<PMDEC,GEOM> C(*this,K); C.check_is_triangulation(); 00401 #endif 00402 } 00403 00404 }; // Constrained_triang_traits<PMDEC,GEOM,NEWEDGE> 00405 00406 CGAL_END_NAMESPACE 00407 00408 #if defined(BOOST_MSVC) 00409 # pragma warning(pop) 00410 #endif 00411 00412 #endif // CGAL_PM_CONSTR_TRIANG_TRAITS_H 00413