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_2/include/CGAL/Nef_2/Constrained_triang_traits.h $
7 // $Id: 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_PM_CONSTR_TRIANG_TRAITS_H
13 #define CGAL_PM_CONSTR_TRIANG_TRAITS_H
14 
15 #include <CGAL/license/Nef_2.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_2/PM_checker.h>
22 #include <cstdlib>
23 #include <string>
24 #include <map>
25 #include <set>
26 #undef CGAL_NEF_DEBUG
27 #define CGAL_NEF_DEBUG 19
28 #include <CGAL/Nef_2/debug.h>
29 
30 #if defined(BOOST_MSVC)
31 #  pragma warning(push)
32 #  pragma warning(disable:4355) // complaint about using 'this' to
33 #endif                          // initialize a member
34 
35 namespace CGAL {
36 
37 struct Do_nothing {
Do_nothingDo_nothing38 Do_nothing() {}
39 template <typename ARG>
operatorDo_nothing40 void operator()(ARG&) const {}
41 };
42 
43 
44 template <typename PMDEC, typename GEOM,
45           typename NEWEDGE = Do_nothing>
46 class Constrained_triang_traits : public PMDEC {
47 public:
48   typedef Constrained_triang_traits<PMDEC,GEOM,NEWEDGE> Self;
49   typedef PMDEC                                         Base;
50 
51   // the types interfacing the sweep:
52   typedef NEWEDGE                   INPUT;
53   typedef typename PMDEC::Plane_map OUTPUT;
54   typedef GEOM                      GEOMETRY;
55 
56   typedef typename GEOM::Point_2     Point;
57   typedef typename GEOM::Segment_2   Segment;
58   typedef typename GEOM::Direction_2 Direction;
59 
60   typedef typename Base::Halfedge_handle   Halfedge_handle;
61   typedef typename Base::Vertex_handle     Vertex_handle;
62   typedef typename Base::Face_handle       Face_handle;
63   typedef typename Base::Halfedge_iterator Halfedge_iterator;
64   typedef typename Base::Vertex_iterator   Vertex_iterator;
65   typedef typename Base::Face_iterator     Face_iterator;
66   typedef typename Base::Halfedge_base     Halfedge_base;
67   typedef typename Base::Halfedge_around_vertex_circulator
68           Halfedge_around_vertex_circulator;
69 
70 
71   using Base::point;
72   using Base::is_isolated;
73   using Base::first_out_edge;
74   using Base::last_out_edge;
75   using Base::source;
76   using Base::target;
77   using Base::twin;
78   using Base::next;
79   using Base::previous;
80   using Base::cyclic_adj_succ;
81   using Base::cyclic_adj_pred;
82   using Base::delete_vertex;
83   using Base::make_first_out_edge;
84 
85   class lt_edges_in_sweepline : public PMDEC
86   {  const Point& p;
87      const Halfedge_handle& e_bottom;
88      const Halfedge_handle& e_top;
89      const GEOMETRY& K;
90     using PMDEC::point;
91     using PMDEC::source;
92     using PMDEC::target;
93   public:
lt_edges_in_sweepline(const Point & pi,const Halfedge_handle & e1,const Halfedge_handle & e2,const PMDEC & D,const GEOMETRY & k)94   lt_edges_in_sweepline(const Point& pi,
95      const Halfedge_handle& e1, const Halfedge_handle& e2,
96      const PMDEC& D, const GEOMETRY& k) :
97        PMDEC(D), p(pi), e_bottom(e1), e_top(e2), K(k) {}
98 
lt_edges_in_sweepline(const lt_edges_in_sweepline & lt)99   lt_edges_in_sweepline(const lt_edges_in_sweepline& lt) :
100      PMDEC(lt), p(lt.p), e_bottom(lt.e_bottom), e_top(lt.e_top), K(lt.K) {}
101 
seg(const Halfedge_handle & e)102   Segment seg(const Halfedge_handle& e) const
103   { return K.construct_segment(point(source(e)),point(target(e))); }
104 
orientation(Halfedge_handle e,const Point & p)105   int orientation(Halfedge_handle e, const Point& p) const
106   { return K.orientation(point(source(e)),point(target(e)),p); }
107 
operator()108   bool operator()(const Halfedge_handle& e1, const Halfedge_handle& e2) const
109   { // Precondition:
110     // [[p]] is identical to the source of either [[e1]] or [[e2]].
111     if (e1 == e_bottom || e2 == e_top) return true;
112     if (e2 == e_bottom || e1 == e_top) return false;
113     if ( e1 == e2 ) return 0;
114     int s = 0;
115     if ( p == point(source(e1)) )      s =   orientation(e2,p);
116     else if ( p == point(source(e2)) ) s = - orientation(e1,p);
117     else CGAL_error_msg("compare error in sweep.");
118     if ( s || source(e1) == target(e1) || source(e2) == target(e2) )
119       return ( s < 0 );
120     s = orientation(e2,point(target(e1)));
121     if (s==0) CGAL_error_msg("parallel edges not allowed.");
122     return ( s < 0 );
123   }
124 
125 
126   }; // lt_edges_in_sweepline
127 
128   class lt_pnts_xy : public PMDEC
129   { const GEOMETRY& K;
130   public:
131     using PMDEC::point;
132 
lt_pnts_xy(const PMDEC & D,const GEOMETRY & k)133    lt_pnts_xy(const PMDEC& D, const GEOMETRY& k) : PMDEC(D), K(k) {}
lt_pnts_xy(const lt_pnts_xy & lt)134    lt_pnts_xy(const lt_pnts_xy& lt) : PMDEC(lt), K(lt.K) {}
operator()135    bool operator()(const Vertex_handle& v1, const Vertex_handle& v2) const
136    { return K.compare_xy(point(v1),point(v2)) < 0; }
137   }; // lt_pnts_xy
138 
139 
140     typedef std::map<Halfedge_handle, Halfedge_handle, lt_edges_in_sweepline>
141             Sweep_status_structure;
142     typedef typename Sweep_status_structure::iterator   ss_iterator;
143     typedef typename Sweep_status_structure::value_type ss_pair;
144     typedef std::set<Vertex_iterator,lt_pnts_xy> Event_Q;
145     typedef typename Event_Q::const_iterator event_iterator;
146 
147     const GEOMETRY&         K;
148     Event_Q                 event_Q;
149     event_iterator          event_it;
150     Vertex_handle           event;
151     Point                   p_sweep;
152     Sweep_status_structure  SL;
153     CGAL::Unique_hash_map<Halfedge_handle,ss_iterator> SLItem;
154     const NEWEDGE&          Treat_new_edge;
155     Halfedge_handle         e_low,e_high; // framing edges !
156     Halfedge_handle         e_search;
157 
Constrained_triang_traits(const INPUT & in,OUTPUT & out,const GEOMETRY & k)158     Constrained_triang_traits(const INPUT& in, OUTPUT& out, const GEOMETRY& k)
159       : Base(out), K(k), event_Q(lt_pnts_xy(*this,K)),
160         SL(lt_edges_in_sweepline(p_sweep,e_low,e_high,*this,K)),
161         SLItem(SL.end()),  Treat_new_edge(in)
162     { CGAL_NEF_TRACEN("Constrained Triangulation Sweep"); }
163 
164 
new_bi_edge(Vertex_handle v1,Vertex_handle v2)165   Halfedge_handle new_bi_edge(Vertex_handle v1, Vertex_handle v2)
166   { // appended at v1 and v2 adj list
167     Halfedge_handle e = Base::new_halfedge_pair(v1,v2);
168     Treat_new_edge(e);
169     return e;
170   }
171 
new_bi_edge(Halfedge_handle e_bf,Halfedge_handle e_af)172   Halfedge_handle new_bi_edge(Halfedge_handle e_bf, Halfedge_handle e_af)
173   { // ccw before e_bf and after e_af
174     Halfedge_handle e = Base::new_halfedge_pair(e_bf,e_af,Halfedge_base(),
175       Base::BEFORE, Base::AFTER);
176     Treat_new_edge(e);
177     return e;
178   }
179 
new_bi_edge(Vertex_handle v,Halfedge_handle e_bf)180   Halfedge_handle new_bi_edge(Vertex_handle v, Halfedge_handle e_bf)
181   { // appended at v's adj list and before e_bf
182     Halfedge_handle e = Base::new_halfedge_pair(v,e_bf,Halfedge_base(),
183       Base::BEFORE);
184     Treat_new_edge(e);
185     return e;
186   }
187 
seg(Halfedge_handle e)188   Segment seg(Halfedge_handle e) const
189   { return K.construct_segment(point(source(e)),point(target(e))); }
190 
dir(Halfedge_handle e)191   Direction dir(Halfedge_handle e) const
192   { return K.construct_direction(point(source(e)),point(target(e))); }
193 
is_forward(Halfedge_handle e)194   bool is_forward(Halfedge_handle e) const
195   { return K.compare_xy(point(source(e)),point(target(e))) < 0; }
196 
197 
198 
199 
edge_is_visible_from(Vertex_handle v,Halfedge_handle e)200   bool edge_is_visible_from(Vertex_handle v, Halfedge_handle e)
201   {
202     Point p =  point(v);
203     Point p1 = point(source(e));
204     Point p2 = point(target(e));
205     return ( K.orientation(p1,p2,p)>0 ); // left_turn
206   }
207 
triangulate_up(Halfedge_handle & e_apex)208   void triangulate_up(Halfedge_handle& e_apex)
209   {
210     CGAL_NEF_TRACEN("triangulate_up "<<seg(e_apex));
211     Vertex_handle v_apex = source(e_apex);
212     while (true) {
213       Halfedge_handle e_vis = previous(twin(e_apex));
214       bool in_sweep_line = (SLItem[e_vis] != SL.end());
215       bool not_visible = !edge_is_visible_from(v_apex,e_vis);
216         CGAL_NEF_TRACEN(" checking "<<in_sweep_line<<not_visible<<" "<<seg(e_vis));
217       if ( in_sweep_line || not_visible) {
218         CGAL_NEF_TRACEN("  STOP"); return;
219       }
220       Halfedge_handle e_back = new_bi_edge(e_apex,e_vis);
221       if ( !is_forward(e_vis) ) make_first_out_edge(twin(e_back));
222       e_apex = e_back;
223       CGAL_NEF_TRACEN(" produced " << seg(e_apex));
224     }
225   }
226 
triangulate_down(Halfedge_handle & e_apex)227   void triangulate_down(Halfedge_handle& e_apex)
228   {
229     CGAL_NEF_TRACEN("triangulate_down "<<seg(e_apex));
230     Vertex_handle v_apex = source(e_apex);
231     while (true) {
232       Halfedge_handle e_vis = next(e_apex);
233       bool in_sweep_line = (SLItem[e_vis] != SL.end());
234       bool not_visible = !edge_is_visible_from(v_apex,e_vis);
235         CGAL_NEF_TRACEN(" checking "<<in_sweep_line<<not_visible<<" "<<seg(e_vis));
236       if ( in_sweep_line || not_visible) {
237           CGAL_NEF_TRACEN("  STOP"); return;
238       }
239       Halfedge_handle e_vis_rev = twin(e_vis);
240       Halfedge_handle e_forw = new_bi_edge(e_vis_rev,e_apex);
241       e_apex = twin(e_forw);
242       CGAL_NEF_TRACEN(" produced " << seg(e_apex));
243     }
244   }
245 
triangulate_between(Halfedge_handle e_upper,Halfedge_handle e_lower)246   void triangulate_between(Halfedge_handle e_upper, Halfedge_handle e_lower)
247   {
248     // we triangulate the interior of the whole chain between
249     // target(e_upper) and target(e_lower)
250     CGAL_assertion(source(e_upper)==source(e_lower));
251     CGAL_NEF_TRACE("triangulate_between\n   "<<seg(e_upper));
252     CGAL_NEF_TRACEN("\n   "<<seg(e_lower));
253     Halfedge_handle e_end = twin(e_lower);
254     while (true) {
255       Halfedge_handle e_vis =  next(e_upper);
256       Halfedge_handle en_vis = next(e_vis);
257       CGAL_NEF_TRACEN(" working on base e_vis " << seg(e_vis));
258       CGAL_NEF_TRACEN(" next is " << seg(en_vis));
259       if (en_vis == e_end) return;
260       e_upper = twin(new_bi_edge(twin(e_vis),e_upper));
261       CGAL_NEF_TRACEN(" produced " << seg(e_upper));
262     }
263   }
264 
process_event()265   void process_event()
266   {
267       CGAL_NEF_TRACEN("\nPROCESS_EVENT " << p_sweep);
268     Halfedge_handle e, ep, eb_low, eb_high, e_end;
269     if ( !is_isolated(event) ) {
270       e = last_out_edge(event);
271       ep = first_out_edge(event);
272     }
273     ss_iterator sit_pred, sit;
274     /* PRECONDITION:
275        only ingoing => e is lowest in ingoing bundle
276        only outgoing => e is highest in outgoing bundle
277        ingoing and outgoing => e is lowest in ingoing bundle */
278     eb_high = e_end = ep;
279     eb_low = e;
280     CGAL_NEF_TRACEN("determining handle in SL");
281     if ( e != Halfedge_handle() ) {
282       point(target(e_search)) = p_sweep; // degenerate loop edge
283       sit_pred = SLItem[e];
284       if ( sit_pred != SL.end())  sit = --sit_pred;
285       else  sit = sit_pred = --SL.upper_bound(e_search);
286     } else { // event is isolated vertex
287       point(target(e_search)) = p_sweep; // degenerate loop edge
288       sit_pred = --SL.upper_bound(e_search);
289     }
290 
291     bool ending_edges(0), starting_edges(0);
292     while ( e != Halfedge_handle() ) { // walk adjacency list clockwise
293       if ( SLItem[e] != SL.end() )
294       {
295         CGAL_NEF_TRACEN("ending " << seg(e));
296         if (ending_edges) triangulate_between(e,cyclic_adj_succ(e));
297         ending_edges = true;
298         SL.erase(SLItem[e]);
299         link_bi_edge_to(e,SL.end());
300         // not in SL anymore
301       }
302 
303       else
304       {
305         CGAL_NEF_TRACEN("starting "<<seg(e));
306         sit = SL.insert(sit,ss_pair(e,e));
307         link_bi_edge_to(e,sit);
308         if ( !starting_edges ) eb_high = cyclic_adj_succ(e);
309         starting_edges = true;
310       }
311 
312       if (e == e_end) break;
313       e = cyclic_adj_pred(e);
314     }
315     if (!ending_edges)
316     {
317       Halfedge_handle e_vis = sit_pred->second;
318       Halfedge_handle e_vis_n = cyclic_adj_succ(e_vis);
319       eb_low = eb_high = new_bi_edge(event,e_vis_n);
320       CGAL_NEF_TRACEN(" producing link "<<seg(eb_low)<<"\n    before "<<seg(e_vis_n));
321     }
322 
323 
324 
325 
326     triangulate_up(eb_high);
327     triangulate_down(eb_low);
328     sit_pred->second = eb_low;
329   }
330 
event_exists()331   bool event_exists()
332   { if ( event_it != event_Q.end() ) {
333       // event is set at end of loop and in init
334       event = *event_it;
335       p_sweep = point(event);
336       return true;
337     }
338     return false;
339   }
340 
procede_to_next_event()341   void procede_to_next_event()
342   { ++event_it; }
343 
link_bi_edge_to(Halfedge_handle e,ss_iterator sit)344   void link_bi_edge_to(Halfedge_handle e, ss_iterator sit) {
345     SLItem[e] = SLItem[twin(e)] = sit;
346   }
347 
initialize_structures()348   void initialize_structures()
349   {
350       CGAL_NEF_TRACEN("initialize_structures ");
351 
352     for ( event=this->vertices_begin(); event != this->vertices_end(); ++event )
353       event_Q.insert(event); // sorted order of vertices
354 
355     event_it = event_Q.begin();
356     if ( event_Q.empty() ) return;
357     event = *event_it;
358     p_sweep = point(event);
359     if ( !is_isolated(event) ) {
360       Halfedge_around_vertex_circulator
361         e(first_out_edge(event)), eend(e);
362       CGAL_For_all(e,eend) {
363         CGAL_NEF_TRACEN("init with "<<PE(e));
364         ss_iterator sit = SL.insert(ss_pair(e,e)).first;
365         link_bi_edge_to(e,sit);
366       }
367     }
368 
369 
370     Vertex_handle v_tmp = this->new_vertex(); point(v_tmp) = Point();
371     e_high = Base::new_halfedge_pair(event,v_tmp);
372     e_low  = Base::new_halfedge_pair(event,v_tmp);
373     // this are two symbolic edges just accessed as sentinels
374     // they carry no geometric information
375     e_search = Base::new_halfedge_pair(v_tmp,v_tmp);
376     // this is just a loop used for searches in SL
377 
378     ss_iterator sit_high = SL.insert(ss_pair(e_high,e_high)).first;
379     ss_iterator sit_low  = SL.insert(ss_pair(e_low,e_low)).first;
380     // inserting sentinels into SL
381     link_bi_edge_to(e_high, sit_high);
382     link_bi_edge_to(e_low , sit_low);
383     // we mark them being in the sweepline, which they will never leave
384 
385 
386     // we move to the second vertex:
387     procede_to_next_event();
388     event_exists(); // sets p_sweep for check invariants
389     CGAL_NEF_TRACEN("EOF initialization");
390   }
391 
complete_structures()392   void complete_structures()
393   {
394     if (e_low != Halfedge_handle()) {
395       delete_vertex(target(e_search));
396     } // removing sentinels and e_search
397   }
398 
399 
check_ccw_local_embedding()400   void check_ccw_local_embedding() const
401   { PM_checker<PMDEC,GEOM> C(*this,K);
402     C.check_order_preserving_embedding(event);
403   }
404 
check_invariants()405   void check_invariants()
406   {
407   #ifdef CGAL_CHECK_EXPENSIVE
408     if ( event_it == event_Q.end() ) return;
409     check_ccw_local_embedding();
410   #endif
411   }
412 
check_final()413   void check_final()
414   {
415   #ifdef CGAL_CHECK_EXPENSIVE
416     PM_checker<PMDEC,GEOM> C(*this,K); C.check_is_triangulation();
417   #endif
418   }
419 
420 }; // Constrained_triang_traits<PMDEC,GEOM,NEWEDGE>
421 
422 } //namespace CGAL
423 
424 #if defined(BOOST_MSVC)
425 #  pragma warning(pop)
426 #endif
427 
428 #endif // CGAL_PM_CONSTR_TRIANG_TRAITS_H
429