1 // Copyright (c) 1997-2000  Max-Planck-Institute Saarbruecken (Germany).
2 // All rights reserved.
3 //
4 // This file is part of CGAL (www.cgal.org).
5 //
6 // $URL: https://github.com/CGAL/cgal/blob/v5.3/Nef_S2/include/CGAL/Nef_S2/SM_constrained_triang_traits.h $
7 // $Id: SM_constrained_triang_traits.h 0779373 2020-03-26T13:31:46+01:00 Sébastien Loriot
8 // SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
9 //
10 //
11 // Author(s)     : Michael Seel <seel@mpi-sb.mpg.de>
12 #ifndef CGAL_SM_CONSTRAINED_TRIANG_TRAITS_H
13 #define CGAL_SM_CONSTRAINED_TRIANG_TRAITS_H
14 
15 #include <CGAL/license/Nef_S2.h>
16 
17 
18 #include <CGAL/basic.h>
19 #include <CGAL/Unique_hash_map.h>
20 #include <CGAL/generic_sweep.h>
21 //#include <CGAL/Nef_S2/SM_checker.h>
22 #include <string>
23 #include <map>
24 #include <set>
25 #undef CGAL_NEF_DEBUG
26 #define CGAL_NEF_DEBUG 139
27 #include <CGAL/Nef_2/debug.h>
28 
29 namespace CGAL {
30 
31 /* For a detailed documentation see the MPI research report 2001-1-003
32    which documents the planar flavor of this baby; only minor deviations
33    are included in this code */
34 
35 template <typename Decorator_, typename Kernel_>
36 class SM_constrained_triang_traits : public Decorator_ {
37 public:
38   typedef SM_constrained_triang_traits<Decorator_,Kernel_> Self;
39   typedef Decorator_                                       Base;
40 
41 
42   typedef typename Kernel_::Point_2     Point;
43   typedef typename Kernel_::Segment_2   Segment;
44 
45   typedef typename Base::SHalfedge_handle   SHalfedge_handle;
46   typedef typename Base::SVertex_handle     SVertex_handle;
47   typedef typename Base::SFace_handle       SFace_handle;
48   typedef typename Base::SHalfedge_iterator SHalfedge_iterator;
49   typedef typename Base::SVertex_iterator   SVertex_iterator;
50   typedef typename Base::SFace_iterator     SFace_iterator;
51   typedef typename Base::SHalfedge_around_svertex_circulator
52           SHalfedge_around_svertex_circulator;
53 
54   // the types interfacing the sweep:
55   typedef std::pair<SVertex_iterator,SVertex_iterator> INPUT;
56   typedef Decorator_                                 OUTPUT;
57   typedef Kernel_                                    GEOMETRY;
58 
59   using Base::is_isolated;
60   using Base::first_out_edge;
61   using Base::last_out_edge;
62   using Base::new_svertex;
63   using Base::assoc_info;
64   using Base::incident_mark;
65   using Base::cyclic_adj_succ;
66   using Base::cyclic_adj_pred;
67   using Base::delete_vertex;
68 
69   class lt_edges_in_sweepline : public Decorator_
70   {  const Point& p;
71      const SHalfedge_handle& e_bottom;
72      const SHalfedge_handle& e_top;
73      const Kernel_& K;
74   public:
lt_edges_in_sweepline(const Point & pi,const SHalfedge_handle & e1,const SHalfedge_handle & e2,const Decorator_ & D,const Kernel_ & k)75   lt_edges_in_sweepline(const Point& pi,
76      const SHalfedge_handle& e1, const SHalfedge_handle& e2,
77      const Decorator_& D, const Kernel_& k) :
78        Decorator_(D), p(pi), e_bottom(e1), e_top(e2), K(k) {}
79 
lt_edges_in_sweepline(const lt_edges_in_sweepline & lt)80   lt_edges_in_sweepline(const lt_edges_in_sweepline& lt) :
81      Decorator_(lt), p(lt.p),
82      e_bottom(lt.e_bottom), e_top(lt.e_top), K(lt.K) {}
83 
seg(const SHalfedge_handle & e)84   Segment seg(const SHalfedge_handle& e) const
85   { return K.construct_segment(point(source(e)),point(target(e))); }
86 
orientation(SHalfedge_handle e,const Point & p)87   int orientation(SHalfedge_handle e, const Point& p) const
88   { return K.orientation(e->source()->point(),e->twin()->source()->point(),p); }
89 
operator()90   bool operator()(const SHalfedge_handle& e1, const SHalfedge_handle& e2) const
91   { // Precondition:
92     // [[p]] is identical to the source of either [[e1]] or [[e2]].
93     if (e1 == e_bottom || e2 == e_top) return true;
94     if (e2 == e_bottom || e1 == e_top) return false;
95     if ( e1 == e2 ) return 0;
96     int s = 0;
97     if ( p == e1->source()->point() )
98       s =   orientation(e2,p);
99     else if ( p == e2->source()->point() )
100       s = - orientation(e1,p);
101     else CGAL_error_msg("compare error in sweep.");
102     if ( s || e1->source() == e1->twin()->source() ||
103          e2->source() == e2->twin()->source())
104       return ( s < 0 );
105     s = orientation(e2,e1->twin()->source()->point());
106     if (s==0) CGAL_error_msg("parallel edges not allowed.");
107     return ( s < 0 );
108   }
109 
110 
111   }; // lt_edges_in_sweepline
112 
113   class lt_pnts_xy : public Decorator_
114   { const Kernel_& K;
115   public:
lt_pnts_xy(const Decorator_ & D,const Kernel_ & k)116    lt_pnts_xy(const Decorator_& D, const Kernel_& k) : Decorator_(D), K(k) {}
lt_pnts_xy(const lt_pnts_xy & lt)117    lt_pnts_xy(const lt_pnts_xy& lt) : Decorator_(lt), K(lt.K) {}
operator()118    int operator()(const SVertex_handle& v1, const SVertex_handle& v2) const
119    { return K.compare_xy(v1->point(),v2->point()) < 0; }
120   }; // lt_pnts_xy
121 
122 
123   typedef std::map<SHalfedge_handle, SHalfedge_handle, lt_edges_in_sweepline>
124           Sweep_status_structure;
125   typedef typename Sweep_status_structure::iterator   ss_iterator;
126   typedef typename Sweep_status_structure::value_type ss_pair;
127   typedef std::set<SVertex_iterator,lt_pnts_xy> Event_Q;
128   typedef typename Event_Q::const_iterator event_iterator;
129 
130   const GEOMETRY&         K;
131   Event_Q                 event_Q;
132   event_iterator          event_it;
133   SVertex_handle           event;
134   Point                   p_sweep;
135   Sweep_status_structure  SL;
136   CGAL::Unique_hash_map<SHalfedge_handle,ss_iterator> SLItem;
137   SHalfedge_handle         e_low,e_high; // framing edges !
138   SHalfedge_handle         e_search;
139   SVertex_iterator         v_first, v_beyond;
140 
SM_constrained_triang_traits(const INPUT & in,OUTPUT & out,const GEOMETRY & k)141   SM_constrained_triang_traits(const INPUT& in, OUTPUT& out,
142                                const GEOMETRY& k)
143       : Base(out), K(k), event_Q(lt_pnts_xy(*this,K)),
144         SL(lt_edges_in_sweepline(p_sweep,e_low,e_high,*this,K)),
145         SLItem(SL.end()), v_first(in.first), v_beyond(in.second)
146     { CGAL_NEF_TRACEN("Constrained Triangulation Sweep"); }
147 
148   /* |treat_new_sedge| is used to forward information that exists
149      at input edges of the triangulation as such it spreads input
150      information to the newly created edges of the triangulation;
151      the used operation incident_mark refers to the base class of
152      |*this| */
treat_new_edge(SHalfedge_handle e)153   void treat_new_edge(SHalfedge_handle e)
154   { assoc_info(e);
155     e->mark() = incident_mark(e) = incident_mark(e->twin()) =
156       incident_mark(e->snext());
157     CGAL_NEF_TRACEN(" treat_new_edge "<<PH(e));
158   }
159 
new_bi_edge(SVertex_handle v1,SVertex_handle v2)160   SHalfedge_handle new_bi_edge(SVertex_handle v1, SVertex_handle v2)
161   { // appended at v1 and v2 adj list
162     SHalfedge_handle e = Base::new_shalfedge_pair(v1,v2);
163     treat_new_edge(e);
164     return e;
165   }
166 
new_bi_edge(SHalfedge_handle e_bf,SHalfedge_handle e_af)167   SHalfedge_handle new_bi_edge(SHalfedge_handle e_bf, SHalfedge_handle e_af)
168   { // ccw before e_bf and after e_af
169     SHalfedge_handle e =
170       Base::new_shalfedge_pair(e_bf,e_af,Base::BEFORE, Base::AFTER);
171     treat_new_edge(e);
172     return e;
173   }
174 
new_bi_edge(SVertex_handle v,SHalfedge_handle e_bf)175   SHalfedge_handle new_bi_edge(SVertex_handle v, SHalfedge_handle e_bf)
176   { // appended at v's adj list and before e_bf
177     SHalfedge_handle e = Base::new_shalfedge_pair(v,e_bf,Base::BEFORE);
178     treat_new_edge(e);
179     return e;
180   }
181 
segment(SHalfedge_handle e)182   Segment segment(SHalfedge_handle e) const
183   { return K.construct_segment(e->source()->point(),e->twin()->source()->point()); }
184 
is_forward(SHalfedge_handle e)185   bool is_forward(SHalfedge_handle e) const
186   { return K.compare_xy(point(source(e)),point(target(e))) < 0; }
187 
edge_is_visible_from(SVertex_handle v,SHalfedge_handle e)188   bool edge_is_visible_from(SVertex_handle v, SHalfedge_handle e)
189   {
190     Point p =  v->point();
191     Point p1 = e->source()->point();
192     Point p2 = e->twin()->source()->point();
193     return ( K.orientation(p1,p2,p) > 0 ); // left_turn
194   }
195 
triangulate_up(SHalfedge_handle & e_apex)196   void triangulate_up(SHalfedge_handle& e_apex)
197   {
198     CGAL_NEF_TRACEN("triangulate_up "<<segment(e_apex));
199     SVertex_handle v_apex = e_apex->source();
200     while (true) {
201       SHalfedge_handle e_vis = e_apex->twin()->sprev();
202       bool in_sweep_line = (SLItem[e_vis] != SL.end());
203       bool not_visible = !edge_is_visible_from(v_apex,e_vis);
204         CGAL_NEF_TRACEN(" checking "<<in_sweep_line<<not_visible<<" "<<segment(e_vis));
205       if ( in_sweep_line || not_visible) {
206         CGAL_NEF_TRACEN("  STOP"); return;
207       }
208       SHalfedge_handle e_back = new_bi_edge(e_apex,e_vis);
209       e_apex = e_back;
210       CGAL_NEF_TRACEN(" produced " << segment(e_apex));
211     }
212   }
213 
triangulate_down(SHalfedge_handle & e_apex)214   void triangulate_down(SHalfedge_handle& e_apex)
215   {
216     CGAL_NEF_TRACEN("triangulate_down "<<segment(e_apex));
217     SVertex_handle v_apex = e_apex->source();
218     while (true) {
219       SHalfedge_handle e_vis = e_apex->snext();
220       bool in_sweep_line = (SLItem[e_vis] != SL.end());
221       bool not_visible = !edge_is_visible_from(v_apex,e_vis);
222         CGAL_NEF_TRACEN(" checking "<<in_sweep_line<<not_visible<<" "<<segment(e_vis));
223       if ( in_sweep_line || not_visible) {
224           CGAL_NEF_TRACEN("  STOP"); return;
225       }
226       SHalfedge_handle e_vis_rev = e_vis->twin();
227       SHalfedge_handle e_forw = new_bi_edge(e_vis_rev,e_apex);
228       e_apex = e_forw->twin();
229       CGAL_NEF_TRACEN(" produced " << segment(e_apex));
230     }
231   }
232 
triangulate_between(SHalfedge_handle e_upper,SHalfedge_handle e_lower)233   void triangulate_between(SHalfedge_handle e_upper, SHalfedge_handle e_lower)
234   {
235     // we triangulate the interior of the whole chain between
236     // target(e_upper) and target(e_lower)
237     CGAL_assertion(e_upper->source()==e_lower->source());
238     CGAL_NEF_TRACE("triangulate_between\n   "<<segment(e_upper));
239     CGAL_NEF_TRACEN("\n   "<<segment(e_lower));
240     SHalfedge_handle e_end = e_lower->twin();
241     while (true) {
242       SHalfedge_handle e_vis =  e_upper->snext();
243       SHalfedge_handle en_vis = e_vis->snext();
244       CGAL_NEF_TRACEN(" working on base e_vis " << segment(e_vis));
245       CGAL_NEF_TRACEN(" next is " << segment(en_vis));
246       if (en_vis == e_end) return;
247       e_upper = new_bi_edge(e_vis->twin(),e_upper)->twin();
248       CGAL_NEF_TRACEN(" produced " << segment(e_upper));
249     }
250   }
251 
process_event()252   void process_event()
253   {
254       CGAL_NEF_TRACEN("\nPROCESS_EVENT " << p_sweep);
255     SHalfedge_handle e, ep, eb_low, eb_high, e_end;
256     if ( !is_isolated(event) ) {
257       e = last_out_edge(event);
258       ep = first_out_edge(event);
259     }
260     ss_iterator sit_pred, sit;
261     /* PRECONDITION:
262        only ingoing => e is lowest in ingoing bundle
263        only outgoing => e is highest in outgoing bundle
264        ingoing and outgoing => e is lowest in ingoing bundle */
265     eb_high = e_end = ep;
266     eb_low = e;
267     CGAL_NEF_TRACEN("determining handle in SL");
268     if ( e != SHalfedge_handle() ) {
269       e_search->twin()->source()->point() = p_sweep; // degenerate loop edge
270       sit_pred = SLItem[e];
271       if ( sit_pred != SL.end())  sit = --sit_pred;
272       else  sit = sit_pred = --SL.upper_bound(e_search);
273     } else { // event is isolated vertex
274       e_search->twin()->source()->point() = p_sweep; // degenerate loop edge
275       sit_pred = --SL.upper_bound(e_search);
276     }
277 
278     bool ending_edges(0), starting_edges(0);
279     while ( e != SHalfedge_handle() ) { // walk adjacency list clockwise
280       if ( SLItem[e] != SL.end() )
281       {
282         CGAL_NEF_TRACEN("ending " << segment(e));
283         if (ending_edges) triangulate_between(e,cyclic_adj_succ(e));
284         ending_edges = true;
285         SL.erase(SLItem[e]);
286         link_bi_edge_to(e,SL.end());
287         // not in SL anymore
288       }
289       else
290       {
291         CGAL_NEF_TRACEN("starting "<<segment(e));
292         sit = SL.insert(sit,ss_pair(e,e));
293         link_bi_edge_to(e,sit);
294         if ( !starting_edges ) eb_high = cyclic_adj_succ(e);
295         starting_edges = true;
296       }
297 
298       if (e == e_end) break;
299       e = cyclic_adj_pred(e);
300     }
301     if (!ending_edges)
302     {
303       SHalfedge_handle e_vis = sit_pred->second;
304       SHalfedge_handle e_vis_n = cyclic_adj_succ(e_vis);
305       eb_low = eb_high = new_bi_edge(event,e_vis_n);
306       CGAL_NEF_TRACEN(" producing link "<<segment(eb_low)<<
307              "\n    before "<<segment(e_vis_n));
308     }
309 
310     triangulate_up(eb_high);
311     triangulate_down(eb_low);
312     sit_pred->second = eb_low;
313   }
314 
event_exists()315   bool event_exists()
316   { if ( event_it != event_Q.end() ) {
317       // event is set at end of loop and in init
318       event = *event_it;
319       p_sweep = event->point();
320       return true;
321     }
322     return false;
323   }
324 
procede_to_next_event()325   void procede_to_next_event()
326   { ++event_it; }
327 
link_bi_edge_to(SHalfedge_handle e,ss_iterator sit)328   void link_bi_edge_to(SHalfedge_handle e, ss_iterator sit) {
329     SLItem[e] = SLItem[e->twin()] = sit;
330   }
331 
initialize_structures()332   void initialize_structures()
333   {
334       CGAL_NEF_TRACEN("initialize_structures ");
335 
336     for ( event = v_first; event != v_beyond; ++event )
337       event_Q.insert(event); // sorted order of vertices
338 
339     event_it = event_Q.begin();
340     if ( event_Q.empty() ) return;
341     event = *event_it;
342     p_sweep = event->point();
343     if ( !is_isolated(event) ) {
344       SHalfedge_around_svertex_circulator
345         e(first_out_edge(event)), eend(e);
346       CGAL_For_all(e,eend) {
347         CGAL_NEF_TRACEN("init with "<<PH(e));
348         ss_iterator sit = SL.insert(ss_pair(e,e)).first;
349         link_bi_edge_to(e,sit);
350       }
351     }
352 
353 
354     SVertex_handle v_tmp = new_svertex(Point());
355     e_high = Base::new_shalfedge_pair(event,v_tmp);
356     e_low  = Base::new_shalfedge_pair(event,v_tmp);
357     // this are two symbolic edges just accessed as sentinels
358     // they carry no geometric information
359     e_search = Base::new_shalfedge_pair(v_tmp,v_tmp);
360     // this is just a loop used for searches in SL
361 
362     ss_iterator sit_high = SL.insert(ss_pair(e_high,e_high)).first;
363     ss_iterator sit_low  = SL.insert(ss_pair(e_low,e_low)).first;
364     // inserting sentinels into SL
365     link_bi_edge_to(e_high, sit_high);
366     link_bi_edge_to(e_low , sit_low);
367     // we mark them being in the sweepline, which they will never leave
368 
369 
370     // we move to the second vertex:
371     procede_to_next_event();
372     event_exists(); // sets p_sweep for check invariants
373     CGAL_NEF_TRACEN("EOF initialization");
374   }
375 
complete_structures()376   void complete_structures()
377   {
378     if (e_low != SHalfedge_handle()) {
379       delete_vertex(e_search->twin()->source());
380     } // removing sentinels and e_search
381   }
382 
383 
check_ccw_local_embedding()384   void check_ccw_local_embedding() const
385   {
386 #ifdef CGAL_CHECK_EXPENSIVEXXX
387     PM_checker<Decorator_,Kernel_> C(*this,K);
388     C.check_order_preserving_embedding(event);
389 #endif
390   }
391 
check_invariants()392   void check_invariants()
393   {
394 #ifdef CGAL_CHECK_EXPENSIVEXXX
395     if ( event_it == event_Q.end() ) return;
396     check_ccw_local_embedding();
397 #endif
398   }
399 
check_final()400   void check_final()
401   {
402 #ifdef CGAL_CHECK_EXPENSIVEXXX
403     PM_checker<Decorator_,Kernel_> C(*this,K); C.check_is_triangulation();
404 #endif
405   }
406 
407 }; // SM_constrained_triang_traits<Decorator_,Kernel_,New_edge_>
408 
409 
410 } //namespace CGAL
411 #endif // CGAL_SM_CONSTRAINED_TRIANG_TRAITS_H
412