1 // Copyright (c) 2005 Rijksuniversiteit Groningen (Netherlands)
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/Skin_surface_3/include/CGAL/Triangulation_incremental_builder_3.h $
7 // $Id: Triangulation_incremental_builder_3.h 254d60f 2019-10-19T15:23:19+02:00 Sébastien Loriot
8 // SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
9 //
10 //
11 // Author(s)     : Nico Kruithof <Nico@cs.rug.nl>
12 
13 #ifndef CGAL_TDS_INCREMENTAL_BUILDER_3_H
14 #define CGAL_TDS_INCREMENTAL_BUILDER_3_H 1
15 
16 #include <CGAL/license/Skin_surface_3.h>
17 
18 #include <CGAL/basic.h>
19 #include <CGAL/Triangulation_data_structure_3.h>
20 #include <CGAL/array.h>
21 #include <set>
22 #include <list>
23 #include <vector>
24 
25 namespace CGAL {
26 
27 template < class Triangulation_3 >
28 class Triangulation_incremental_builder_3
29 {
30 public:
31   typedef Triangulation_3                    T;
32   typedef typename T::Vertex_handle          Vertex_handle;
33   typedef typename T::Cell_handle            Cell_handle;
34   typedef typename T::Facet                  Facet;
35 
36 public:
37   Triangulation_incremental_builder_3( T &t, bool verbose = false)
t(t)38     : t(t), m_verbose(verbose) {
39   }
40 
~Triangulation_incremental_builder_3()41   ~Triangulation_incremental_builder_3() {}
42 
begin_triangulation(int dim)43   void begin_triangulation(int dim)
44   {
45     t.clear();
46     t.tds().delete_cell(t.infinite_vertex()->cell());
47     // t.infinite = add_vertex();
48     t.tds().set_dimension(dim);
49   }
50 
end_triangulation()51   void end_triangulation()
52   {
53     construct_infinite_cells();
54     CGAL_assertion(t.infinite_vertex()->cell() != Cell_handle());
55   }
56 
57   Vertex_handle add_vertex();
58   Cell_handle add_cell(Vertex_handle vh0, Vertex_handle vh1,
59                        Vertex_handle vh2, Vertex_handle vh3);
60 
61 private:
62   void construct_infinite_cells();
63   Cell_handle add_infinite_cell(Cell_handle ch, int i);
64 
65   void glue_cells(Cell_handle ch0, int ind0, Cell_handle ch1, int ind1);
66 
67   // Interior facets of the simplical cell:
68   typedef std::pair < Vertex_handle, Vertex_handle >    Vpair;
69   typedef std::map < Vpair, Facet >                     MapPair;
70   typedef typename MapPair::iterator                    MapPairIt;
71   typedef std::array < Vertex_handle, 3 >             Vtriple;
72   typedef std::map < Vtriple, Facet >                   MapTriple;
73   typedef typename MapTriple::iterator                  MapTripleIt;
74 
facet(Vertex_handle vh1,Vertex_handle vh2,Vertex_handle vh3)75   Vtriple facet(Vertex_handle vh1, Vertex_handle vh2, Vertex_handle vh3)
76   {
77     if (vh1 < vh2) {
78       if (vh2 < vh3) {
79         return CGAL::make_array(vh1,vh2,vh3);
80       } else {
81         if (vh1 < vh3) {
82           return CGAL::make_array(vh1,vh3,vh2);
83         } else {
84           return CGAL::make_array(vh3,vh1,vh2);
85         }
86       }
87     }
88     if (vh1 < vh3) {
89       return CGAL::make_array(vh2,vh1,vh3);
90     } else {
91       if (vh2 < vh3) {
92         return CGAL::make_array(vh2,vh3,vh1);
93       } else {
94         return CGAL::make_array(vh3,vh2,vh1);
95       }
96     }
97   }
98 
99   MapTriple facets;
100 
101   T &t;
102   bool m_verbose;
103 };
104 
105 template < class TDS_>
106 typename Triangulation_incremental_builder_3< TDS_ >::Vertex_handle
add_vertex()107 Triangulation_incremental_builder_3< TDS_ >::add_vertex() {
108   return t.tds().create_vertex();
109 }
110 
111 template < class TDS_>
112 typename Triangulation_incremental_builder_3< TDS_ >::Cell_handle
add_cell(Vertex_handle vh0,Vertex_handle vh1,Vertex_handle vh2,Vertex_handle vh3)113 Triangulation_incremental_builder_3< TDS_ >::add_cell(
114   Vertex_handle vh0, Vertex_handle vh1, Vertex_handle vh2, Vertex_handle vh3)
115 {
116   CGAL_assertion(vh0 != nullptr); CGAL_assertion(vh1 != nullptr);
117   CGAL_assertion(vh2 != nullptr); CGAL_assertion(vh3 != nullptr);
118   CGAL_assertion(vh0 != vh1); CGAL_assertion(vh0 != vh2); CGAL_assertion(vh0 != vh3);
119   CGAL_assertion(vh1 != vh2); CGAL_assertion(vh1 != vh3); CGAL_assertion(vh2 != vh3);
120 
121   Cell_handle ch =  t.tds().create_cell(vh0, vh1, vh2, vh3);
122   // Neighbors are by default set to nullptr
123   vh0->set_cell(ch); vh1->set_cell(ch);
124   vh2->set_cell(ch); vh3->set_cell(ch);
125 
126   for (int i=0; i<4; i++) {
127     Vtriple vtriple=facet(ch->vertex((i+1)&3),
128                           ch->vertex((i+2)&3),
129                           ch->vertex((i+3)&3));
130 
131     std::pair<MapTripleIt,bool> res = facets.insert(std::make_pair(vtriple, Facet(ch, i)));
132     if (! res.second) { // we found an element with this key
133       MapTripleIt neighbIt = res.first;
134       Facet f = (*neighbIt).second;
135       glue_cells(f.first, f.second, ch, i);
136       facets.erase(neighbIt);
137       CGAL_assertion(f.first->neighbor(f.second) != nullptr);
138       CGAL_assertion(ch->neighbor(i) != nullptr);
139     } else {
140       CGAL_assertion(ch->neighbor(i) == nullptr);
141     }
142   }
143 
144   return ch;
145 }
146 
147 template < class TDS_>
148 typename Triangulation_incremental_builder_3< TDS_ >::Cell_handle
add_infinite_cell(Cell_handle ch0,int i)149 Triangulation_incremental_builder_3< TDS_ >::add_infinite_cell(
150   Cell_handle ch0, int i)
151 {
152   CGAL_assertion(ch0->neighbor(i) == nullptr);
153   Vertex_handle vh[4];
154   vh[i] = t.infinite_vertex();
155   vh[(i+1)&3] = ch0->vertex((i+1)&3);
156   vh[(i+2)&3] = ch0->vertex((i+3)&3);
157   vh[(i+3)&3] = ch0->vertex((i+2)&3);
158   Cell_handle ch1 =  t.tds().create_cell(vh[0], vh[1], vh[2], vh[3]);
159   // Neighbors are set to nullptr
160   // Do not set points to the infinite cell. All finite vertices point to
161   // finite cells.
162   vh[i]->set_cell(ch1);
163   glue_cells(ch0, i, ch1, i);
164   return ch1;
165 }
166 
167 template < class TDS_>
168 void
glue_cells(Cell_handle ch0,int ind0,Cell_handle ch1,int ind1)169 Triangulation_incremental_builder_3< TDS_ >::glue_cells(
170   Cell_handle ch0, int ind0, Cell_handle ch1, int ind1)
171 {
172   CGAL_assertion(ch0->has_vertex(ch1->vertex((ind1+1)&3)));
173   CGAL_assertion(ch0->has_vertex(ch1->vertex((ind1+2)&3)));
174   CGAL_assertion(ch0->has_vertex(ch1->vertex((ind1+3)&3)));
175 
176   CGAL_assertion(ch1->has_vertex(ch0->vertex((ind0+1)&3)));
177   CGAL_assertion(ch1->has_vertex(ch0->vertex((ind0+2)&3)));
178   CGAL_assertion(ch1->has_vertex(ch0->vertex((ind0+3)&3)));
179 
180   CGAL_assertion(ch0->index(ch1->vertex((ind1+1)&3)) != ind0);
181   CGAL_assertion(ch0->index(ch1->vertex((ind1+2)&3)) != ind0);
182   CGAL_assertion(ch0->index(ch1->vertex((ind1+3)&3)) != ind0);
183 
184   CGAL_assertion(ch1->index(ch0->vertex((ind0+1)&3)) != ind1);
185   CGAL_assertion(ch1->index(ch0->vertex((ind0+2)&3)) != ind1);
186   CGAL_assertion(ch1->index(ch0->vertex((ind0+3)&3)) != ind1);
187 
188   ch0->set_neighbor(ind0, ch1);
189   ch1->set_neighbor(ind1, ch0);
190 }
191 
192 // Adds infinite cells to the facets on the convex hull
193 template < class TDS_>
194 void
construct_infinite_cells()195 Triangulation_incremental_builder_3< TDS_ >::construct_infinite_cells()
196 {
197   MapTripleIt ch_facet_it;
198   MapPair     ch_edges;
199   MapPairIt   ch_edge_it;
200   Vertex_handle vh1, vh2;
201 
202   for (ch_facet_it = facets.begin();
203        ch_facet_it != facets.end();
204        ch_facet_it ++) {
205     Facet ch_facet = (*ch_facet_it).second;
206     Cell_handle ch0 = ch_facet.first;
207     int ind0 = ch_facet.second;
208     Cell_handle ch1 = add_infinite_cell(ch0, ind0);
209     // Index of ch1 is also ind0
210     CGAL_assertion(ch0->neighbor(ind0) != nullptr);
211     CGAL_assertion(ch1->neighbor(ind0) != nullptr);
212 
213     for (int i=1; i<4; i++) {
214       int i1 = (i==1?2:1);
215       int i2 = (i==3?2:3);
216       if (ch1->vertex((ind0+i1)&3) < ch1->vertex((ind0+i2)&3)) {
217         vh1 = ch1->vertex((ind0+i1)&3);
218         vh2 = ch1->vertex((ind0+i2)&3);
219       } else {
220         vh1 = ch1->vertex((ind0+i2)&3);
221         vh2 = ch1->vertex((ind0+i1)&3);
222       }
223       ch_edge_it = ch_edges.find(Vpair(vh1,vh2));
224       if (ch_edge_it != ch_edges.end()) {
225         Facet f_opp = (*ch_edge_it).second;
226         glue_cells(f_opp.first, f_opp.second, ch1, (ind0+i)&3);
227         ch_edges.erase(ch_edge_it);
228         CGAL_assertion(f_opp.first->neighbor(f_opp.second) != nullptr);
229         CGAL_assertion(ch1->neighbor((ind0+i)&3) != nullptr);
230       } else {
231         ch_edges[Vpair(vh1,vh2)] = Facet(ch1, (ind0+i)&3);
232         CGAL_assertion(ch1->neighbor((ind0+i)&3) == nullptr);
233       }
234     }
235   }
236   CGAL_assertion(ch_edges.empty());
237 }
238 
239 } //namespace CGAL
240 
241 #endif // CGAL_TDS_INCREMENTAL_BUILDER_3_H //
242