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/Convex_hull_d/include/CGAL/Regular_complex_d.h $
7 // $Id: Regular_complex_d.h d1a323c 2020-03-26T19:24:14+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 //---------------------------------------------------------------------
13 // file generated by notangle from regl_complex.lw
14 // please debug or modify LEDA web file
15 // mails and bugs: Michael.Seel@mpi-sb.mpg.de
16 // based on LEDA architecture by S. Naeher, C. Uhrig
17 // coding: K. Mehlhorn, M. Seel
18 // debugging and templatization: M. Seel
19 //---------------------------------------------------------------------
20 
21 #ifndef CGAL_REGULAR_COMPLEX_D_H
22 #define CGAL_REGULAR_COMPLEX_D_H
23 
24 #include <CGAL/license/Convex_hull_d.h>
25 
26 #define CGAL_DEPRECATED_HEADER "<CGAL/Regular_complex_d.h>"
27 #define CGAL_DEPRECATED_MESSAGE_DETAILS \
28   "The Triangulation package (see https://doc.cgal.org/latest/Triangulation) should be used instead."
29 #include <CGAL/internal/deprecation_warning.h>
30 
31 #include <CGAL/basic.h>
32 #include <CGAL/Iterator_project.h>
33 #include <CGAL/Compact_container.h>
34 #include <vector>
35 #include <list>
36 #include <cstddef>
37 
38 #include <CGAL/Kernel_d/debug.h>
39 
40 #ifdef CGAL_USE_LEDA
41 #include <LEDA/system/memory.h>
42 #endif
43 
44 namespace CGAL {
45 
46 template <class R> class RC_simplex_d;
47 template <class R> class RC_vertex_d;
48 template <class R> class Regular_complex_d;
49 template <class R> class Convex_hull_d;
50 
51 #define forall_rc_vertices(x,RC)\
52 for(x = (RC).vertices_begin(); x != (RC).vertices_end(); ++x)
53 #define forall_rc_simplices(x,RC)\
54 for(x = (RC).simplices_begin(); x != (RC).simplices_end(); ++x)
55 
56 
57 template <class Refs>
58 class RC_vertex_d
59 
60 { typedef RC_vertex_d<Refs> Self;
61   typedef typename Refs::Point_d Point_d;
62   typedef typename Refs::Vertex_handle Vertex_handle;
63   typedef typename Refs::Simplex_handle Simplex_handle;
64   typedef typename Refs::R R;
65   friend class Regular_complex_d<R>;
66   friend class Convex_hull_d<R>;
67   friend class RC_simplex_d<Refs>;
68   Simplex_handle s_;
69   int            index_;
70   Point_d        point_;
71 
set_simplex(Simplex_handle s)72   void set_simplex(Simplex_handle s) { s_=s; }
set_index(int i)73   void set_index(int i) { index_=i; }
set_point(const Point_d & p)74   void set_point(const Point_d& p) { point_=p; }
75 
76 public:
RC_vertex_d(Simplex_handle s,int i,const Point_d & p)77   RC_vertex_d(Simplex_handle s, int i, const Point_d& p) :
78     s_(s), index_(i), point_(p) {}
RC_vertex_d(const Point_d & p)79   RC_vertex_d(const Point_d& p) : point_(p), pp(nullptr) {}
RC_vertex_d()80   RC_vertex_d() :  s_(), pp(nullptr) {}
81   // beware that ass_point was initialized here by nil_point
~RC_vertex_d()82   ~RC_vertex_d() {}
83 
simplex()84   Simplex_handle simplex() const { return s_; }
index()85   int index() const { return index_; }
point()86   const Point_d& point() const { return point_; }
87 
88   void* pp;
for_compact_container()89   void*   for_compact_container() const { return pp; }
for_compact_container(void * p)90   void for_compact_container(void *p) { pp = p; }
91 
92 #ifdef CGAL_USE_LEDA
93   LEDA_MEMORY(RC_vertex_d)
94 #endif
95 };
96 
97 
98 template <class Refs>
99 class RC_simplex_d
100 { typedef RC_simplex_d<Refs>  Self;
101   typedef typename Refs::Point_d Point_d;
102   typedef typename Refs::Vertex_handle Vertex_handle;
103   typedef typename Refs::Simplex_handle Simplex_handle;
104   typedef typename Refs::R R;
105   friend class Regular_complex_d<R>;
106   friend class Convex_hull_d<R>;
107 protected:
108   std::vector<Vertex_handle>   vertices;    // array of vertices
109   std::vector<Simplex_handle>  neighbors;   // opposite simplices
110   std::vector<int>             opposite_vertices;
111                                // indices of opposite vertices
112 
113   //------ only for convex hulls ------------------
114   typedef typename R::Hyperplane_d Hyperplane_d;
115   Hyperplane_d h_base;   // hyperplane supporting base facet
116   bool         visited_; // visited mark when traversing
117   //------ only for convex hulls ------------------
118 
vertex(int i)119   Vertex_handle  vertex(int i) const { return vertices[i]; }
neighbor(int i)120   Simplex_handle neighbor(int i) const { return neighbors[i]; }
opposite_vertex_index(int i)121   int opposite_vertex_index(int i) const { return opposite_vertices[i]; }
122 
set_vertex(int i,Vertex_handle v)123   void set_vertex(int i, Vertex_handle v) { vertices[i] = v; }
set_neighbor(int i,Simplex_handle s)124   void set_neighbor(int i, Simplex_handle s) { neighbors[i]=s; }
set_opposite_vertex_index(int i,int index)125   void set_opposite_vertex_index(int i, int index)
126   { opposite_vertices[i]=index; }
127 
128   //------ only for convex hulls ------------------
hyperplane_of_base_facet()129   Hyperplane_d hyperplane_of_base_facet() const { return h_base; }
set_hyperplane_of_base_facet(const Hyperplane_d & h)130   void set_hyperplane_of_base_facet(const Hyperplane_d& h) { h_base = h; }
visited()131   bool& visited() { return visited_; }
132   //------ only for convex hulls ------------------
133 
134 public:
135   typedef typename std::vector<Vertex_handle>::const_iterator
136           VIV_iterator;
137   struct Point_from_VIV_iterator {
138     typedef Vertex_handle argument_type;
139     typedef Point_d      result_type;
operatorPoint_from_VIV_iterator140     result_type& operator()(argument_type& x) const
141     { return x->point(); }
operatorPoint_from_VIV_iterator142     const result_type& operator()(const argument_type& x) const
143     { return x->point(); }
144   };
145 
146   typedef CGAL::Iterator_project<VIV_iterator,Point_from_VIV_iterator,
147     const Point_d&, const Point_d*> Point_const_iterator;
148 
points_begin()149   Point_const_iterator points_begin() const
150   { return Point_const_iterator(vertices.begin()); }
points_end()151   Point_const_iterator points_end() const
152   { return Point_const_iterator(vertices.end()); }
153 
154   void* pp;
for_compact_container()155   void*   for_compact_container() const { return pp; }
for_compact_container(void * p)156   void for_compact_container(void *p) { pp = p; }
157 
158   #if 0
159   struct Point_const_iterator {
160     typedef Point_const_iterator self;
161     typedef std::random_access_iterator_tag iterator_category;
162     typedef const Point_d&                  value_type;
163     typedef std::ptrdiff_t                  difference_type;
164     typedef const Point_d*                  pointer;
165     typedef const Point_d&                  reference;
166 
167     typedef typename std::vector<Vertex_handle>::const_iterator
168       ra_vertex_iterator;
169 
170     Point_const_iterator() : _it() {}
171     Point_const_iterator(ra_vertex_iterator it) : _it(it) {}
172 
173     value_type operator*() const { return (*_it)->point(); }
174     pointer operator->() const { return &(operator*()); }
175 
176     self& operator++() { ++_it; return *this; }
177     self  operator++(int) { self tmp = *this; ++_it; return tmp; }
178     self& operator--() { --_it; return *this; }
179     self  operator--(int) { self tmp = *this; --_it; return tmp; }
180 
181     self& operator+=(difference_type i) { _it+=i; return *this; }
182     self& operator-=(difference_type i) { _it-=i; return *this; }
183     self operator+(difference_type i) const
184     { self tmp=*this; tmp+=i; return tmp; }
185     self operator-(difference_type i) const
186     { self tmp=*this; tmp-=i; return tmp; }
187 
188     difference_type operator-(self x) const { return _it-x._it; }
189     reference operator[](difference_type i) { return *(*this + i); }
190 
191     bool operator==(const self& x) const { return _it==x._it; }
192     bool operator!=(const self& x) const { return ! (*this==x); }
193     bool operator<(self x) const { (x - *this) > 0; }
194 
195     private:
196       ra_vertex_iterator _it;
197   }; // Point_const_iterator
198 
199   Point_const_iterator points_begin() const
200   { return Point_const_iterator(vertices.begin()); }
201   Point_const_iterator points_end() const
202   { return Point_const_iterator(vertices.end()); }
203 
204   #endif
205 
206 
RC_simplex_d()207   RC_simplex_d() : pp(nullptr) {}
RC_simplex_d(int dmax)208   RC_simplex_d(int dmax) :
209     vertices(dmax+1), neighbors(dmax+1), opposite_vertices(dmax+1), pp(nullptr)
210   { for (int i = 0; i <= dmax; i++) {
211       neighbors[i] = Simplex_handle();
212       vertices[i] = Vertex_handle();
213       opposite_vertices[i] = -1;
214     }
215     visited_ = false;
216   }
~RC_simplex_d()217   ~RC_simplex_d() {}
218 
219   void print(std::ostream& O=std::cout) const
220   {
221     O << "RC_simplex_d {" ;
222     for(int i=0;i < int(vertices.size());++i) {
223       Vertex_handle v = vertices[i];
224       if ( v != Vertex_handle() ) O << v->point();
225       else O << "(nil)";
226     }
227     O << "}";
228   }
229 
230 #ifdef CGAL_USE_LEDA
231   LEDA_MEMORY(RC_simplex_d)
232 #endif
233 
234 };
235 
236 template <typename R>
237 std::ostream& operator<<(std::ostream& O, const RC_simplex_d<R>& s)
238 { s.print(O); return O; }
239 
240 
241 /*{\Manpage {Regular_complex_d}{R}{Regular Simplicial Complex}{C}}*/
242 /*{\Mdefinition
243 
244 An instance |\Mvar| of type |\Mname| is a regular abstract or concrete
245 simplicial complex. An abstract simplicial complex is a family |\Mvar|
246 of subsets of some set $V$, called the vertex set of the complex,
247 which is closed under the subset relation, i.e., if a set $s$ belongs
248 to the family then all its subsets do. A set $s$ of cardinality $k +
249 1$ is called a $k$-simplex and $k$ is called its dimension.  If $s$ is
250 a subset of $t$ then $s$ is called a subsimplex or face of $t$. A
251 vertex $v$ is called incident to a simplex $s$ if $v$ is an element of
252 $s$.  A simplex is called \emph{maximal} if it is not a face of any
253 simplex in |\Mvar|. Two simplices of dimension $k$ are called
254 neighbors if they share $k-1$ vertices. A complex is connected if its
255 set of maximal simplices forms a connected set under the neighboring
256 relation.  A simplicial complex is called \emph{regular} if all
257 maximal simplices in the complex have the same dimension and if the
258 complex is connected.
259 
260 A concrete simplicial complex is an abstract simplicial complex in
261 which a point in some ambient space is associated with each vertex.
262 We use |dim| to denote the dimension of ambient space.  Simplices are
263 now interpreted geometrically as sets of points in ambient space,
264 namely as the convex hulls of (the points associated with) their
265 vertices. A $0$-simplex is a point, a $1$-simplex is a line segment, a
266 $2$-simplex is a triangle, a $3$-simplex is a tetrahedron,
267 etc.. \emph{The simplices is a concrete simplicial complex must
268 satisfy the additional conditions that the points associated with the
269 vertices of any simplex are affinely independet and that the
270 intersection of any two simplices is a face of both.} We will write
271 simplicial complex instead of concrete simplicial complex in the
272 sequel.
273 
274 All maximal simplices in a regular simplicial complex have the same
275 dimension, which we denote |dcur|.  For each maximal
276 simplex\cgalFootnote{we drop the adjective maximal in the sequel} in
277 |\Mvar| there is an item of type |RC_simplex_d| and for each vertex
278 there is an item of type |rc_vertex|.  Each maximal simplex has |1+dcur|
279 vertices indexed from $0$ to |dcur|. For any simplex $s$ and any
280 index $i$, |C.vertex_of(s,i)| returns the $i$-th vertex of $s$. There
281 may or may not be a simplex $t$ opposite to (the vertex with index)
282 $i$, i.e., a maximal simplex $t$ having
283 \{|C.vertex_of(s,0)|,|C.vertex_of(s,1)|,\ldots,
284 |C.vertex_of(s,dcur)|\} - \{|C.vertex_of(s,i)|\} in its vertex
285 set.  The function |C.opposite(s,i)| returns $t$ if it exists and
286 returns |nil| otherwise. If $t$ exists then $s$ and $t$ share |dcur|
287 vertices, namely all but vertex $i$ of $s$ and vertex
288 |C.opposite_vertex(s,i)| of $t$. Assume that $t = |C.opposite(s,i)|$
289 exists and let |j = C.opposite_vertex(s,i)|. Then |s = C.opposite(t,j)|
290 and |i = C.opposite_vertex(t,j)| and
291 \begin{eqnarray*}
292 \lefteqn{\{|C.vertex_of(s,0)|,|C.vertex_of(s,1)|,\ldots,
293 |C.vertex_of(s,dcur)|\} - \{|C.vertex_of(s,i)|\} =} \\ & &
294 \{|C.vertex_of(t,0)|,|C.vertex_of(t,1)|,\ldots,|C.vertex_of(t,dcur)|\}
295 - \{|C.vertex_of(t,j)|\}.  \end{eqnarray*} In general, a
296 vertex belongs to many simplices. For an |rc_vertex| $v$,
297 the functions |C.simplex(v)| and |C.index(v)| return a pair $(s,i)$
298 such that |v = C.vertex_of(s,i)|.
299 
300 The class |regl_complex| has a static member |nil_point| of type
301 |Point_d|. This point is different (= not indentical) from any user
302 defined point and is the point associated with every vertex of an
303 abstract simplicial complex. It simulates the use of |nil| to denote
304 an undefined object.
305 
306 Regular complexes are designed as the base class for triangulations of
307 convex hulls and Delaunay triangulations in higher dimensional
308 space. We have not used them yet for any other purpose. Regular
309 complexes are built by constructing vertices and simplices, by
310 assigning positions to vertices and vertices to simplices, and by
311 establishing neighbor relations. The update operations do not check
312 whether the data structure built actually encodes a simplicial
313 complex.  The class provides a function |is_valid()| that performs a
314 partial check whether the data structure encodes a simplicial
315 complex. It is not checked whether two simplices intersect without
316 sharing a face.
317 }*/
318 
319 template <class R_>
320 class Regular_complex_d
321 {
322 typedef Regular_complex_d<R_> Self;
323 public:
324 /*{\Mtypes 4}*/
325 typedef R_ R;
326 
327 typedef RC_vertex_d<Self> Vertex;
328 
329   typedef CGAL::Compact_container<Vertex> Vertex_list;
330 
331 typedef typename Vertex_list::iterator       Vertex_handle;
332 /*{\Mtypemember the handle type for vertices of the complex.}*/
333 typedef typename Vertex_list::const_iterator Vertex_const_handle;
334 
335 typedef typename Vertex_list::iterator       Vertex_iterator;
336 /*{\Mtypemember the iterator type for vertices of the complex.}*/
337 typedef typename Vertex_list::const_iterator Vertex_const_iterator;
338 
339 typedef RC_simplex_d<Self> Simplex;
340 typedef CGAL::Compact_container<Simplex>    Simplex_list;
341 
342 typedef typename Simplex_list::iterator       Simplex_handle;
343 /*{\Mtypemember the handle type for simplices of the complex.}*/
344 typedef typename Simplex_list::const_iterator Simplex_const_handle;
345 
346 typedef typename Simplex_list::iterator       Simplex_iterator;
347 /*{\Mtypemember the iterator type for simplices of the complex.}*/
348 typedef typename Simplex_list::const_iterator Simplex_const_iterator;
349 
350 typedef typename R::Point_d Point_d;
351 
352 protected:
353   const R& Kernel_;
354   int dcur; // dimension of the current complex
355   int dmax; // dimension of ambient space
356 
357   Vertex_list  vertices_; // list of all vertices
358   Simplex_list simplices_; // list of all simplices
359 
360 /* the default copy constructor and assignment operator for class
361    regl_complex work incorrectly; it is therefore good practice to
362    either implement them correctly or to make them inaccessible. We do
363    the latter. */
364 
365 private:
366 
367   Regular_complex_d(const Regular_complex_d<R>& );
368   Regular_complex_d& operator=(const Regular_complex_d<R>& );
369 
clean_dynamic_memory()370   void clean_dynamic_memory()
371   {
372     vertices_.clear();
373     simplices_.clear();
374 }
375 
376 public:
377 
378 /*{\Mcreation}*/
379 
380 Regular_complex_d(int d = 2, const R& Kernel = R())
381 /*{\Mcreate creates an instance |\Mvar| of type |\Mtype|. The
382 dimension of the underlying space is $d$ and |\Mvar| is initialized to
383 the empty regular complex.  Thus |dcur| equals $-1$. The traits class
384 |R| specifies the models of all types and the implementations of
385 all geometric primitives used by the regular complex class. The
386 |Kernel| parameter allows you to carry fixed geometric information
387 into the data type. For the default kernel traits |Homogeneous_d|
388 the default construction of |Kernel| is enough.
389 
390 In the following we use further template parameters like the point
391 type |Point_d=R::Point_d|.  At this point, it suffices to say that
392 |Point_d| represents points in $d$-space. The complete specification of
393 the traits class is to be found at the end of this manual page.}*/
Kernel_(Kernel)394  : Kernel_(Kernel) { dmax = d; dcur = -1; }
395 
396 
~Regular_complex_d()397 ~Regular_complex_d() { clean_dynamic_memory(); }
398 
399 /* In the destructor for |Regular_complex_d|, we have to release the
400    storage which was allocated for the simplices and the vertices. */
401 
402 /*{\Mtext The data type |\Mtype| offers neither copy constructor nor
403           assignment operator.}*/
404 
405 /*{\Moperations 3 3}*/
406 /*{\Mtext \headerline{Access Operations}}*/
407 
dimension()408 int dimension() const { return dmax; }
409 /*{\Mop returns the dimension of ambient space}*/
410 
current_dimension()411 int current_dimension() const { return dcur; }
412 /*{\Mop returns the current dimension of the simplices in the
413 complex.}*/
414 
vertex(Simplex_handle s,int i)415 Vertex_handle vertex(Simplex_handle s, int i) const
416 /*{\Mop returns the $i$-th vertex of $s$.\\
417 \precond $0 \leq i \leq |current_dimension|$. }*/
418 { CGAL_assertion(0<=i&&i<=dcur);
419   return s->vertex(i); }
420 
vertex(Simplex_const_handle s,int i)421 Vertex_const_handle vertex(Simplex_const_handle s, int i) const
422 { CGAL_assertion(0<=i&&i<=dcur);
423   return s->vertex(i); }
424 
associated_point(Vertex_handle v)425 Point_d associated_point(Vertex_handle v) const
426 /*{\Mop returns the point associated with vertex |v|.}*/
427 { return v->point(); }
428 
associated_point(Vertex_const_handle v)429 Point_d associated_point(Vertex_const_handle v) const
430 { return v->point(); }
431 
432 
index(Vertex_handle v)433 int index(Vertex_handle v) const
434 /*{\Mop returns the index of $v$ in |C.simplex(v)|.}*/
435 { return v->index(); }
436 
index(Vertex_const_handle v)437 int index(Vertex_const_handle v) const
438 { return v->index(); }
439 
440 
simplex(Vertex_handle v)441 Simplex_handle simplex(Vertex_handle v) const
442 /*{\Mop returns a simplex of which $v$ is a vertex. Note that this
443 simplex is not unique. }*/
444 { return v->simplex(); }
445 
simplex(Vertex_const_handle v)446 Simplex_const_handle simplex(Vertex_const_handle v) const
447 { return v->simplex(); }
448 
associated_point(Simplex_handle s,int i)449 Point_d associated_point(Simplex_handle s, int i) const
450 /*{\Mop same as |C.associated_point(C.vertex(s,i))|.}*/
451 { return associated_point(vertex(s,i)); }
452 
associated_point(Simplex_const_handle s,int i)453 Point_d associated_point(Simplex_const_handle s, int i) const
454 { return associated_point(vertex(s,i)); }
455 
opposite_simplex(Simplex_handle s,int i)456 Simplex_handle opposite_simplex(Simplex_handle s,int i) const
457 /*{\Mop returns the simplex opposite to the $i$-th vertex of $s$
458 (|Simplex_handle()| is there is no such simplex).\\ \precond $0 \leq i \leq |dcur|$. }*/
459 { CGAL_assertion(0<=i&&i<=dcur);
460   return s->neighbor(i); }
461 
opposite_simplex(Simplex_const_handle s,int i)462 Simplex_const_handle opposite_simplex(Simplex_const_handle s,int i) const
463 { CGAL_assertion(0<=i&&i<=dcur);
464   return s->neighbor(i); }
465 
index_of_opposite_vertex(Simplex_handle s,int i)466 int index_of_opposite_vertex(Simplex_handle s, int i) const
467 { CGAL_assertion(0<=i&&i<=dcur);
468   return s->opposite_vertex_index(i); }
469 /*{\Mop returns the index of the vertex opposite to the $i$-th vertex
470 of $s$. \precond $0 \leq i \leq |dcur|$ and there is a
471 simplex opposite to the $i$-th vertex of $s$.}*/
472 
index_of_opposite_vertex(Simplex_const_handle s,int i)473 int index_of_opposite_vertex(Simplex_const_handle s, int i) const
474 { CGAL_assertion(0<=i&&i<=dcur);
475   return s->opposite_vertex_index(i); }
476 
477 
478 
479 /*{\Mtext \headerline{Update Operations}
480 
481 We give operations that allow to update a regular complex. They have
482 to be used with care as they may invalidate the data structure.}*/
483 
484 void clear(int d = 0)
485 /*{\Mop reinitializes |\Mvar| to the empty complex in dimension |dim|.}*/
486 { clean_dynamic_memory();
487   dmax = d; dcur = -1;
488 }
489 
set_current_dimension(int d)490 void set_current_dimension(int d) { dcur = d; }
491 /*{\Mop sets |dcur| to |d|. }*/
492 
new_simplex()493 Simplex_handle new_simplex()
494 /*{\Mop adds a new simplex to |\Mvar| and returns it. The new simplex
495 has no vertices yet.}*/
496 {
497   Simplex s(dmax);
498   Simplex_handle  h = simplices_.insert(s);
499   return h;
500 }
501 
new_vertex()502 Vertex_handle  new_vertex()
503 /*{\Mop adds a new vertex to |\Mvar| and returns it. The new vertex
504         has no associated simplex nor index yet. The associated point
505         is the point |Regular_complex_d::nil_point| which is a static
506         member of class |Regular_complex_d.|}*/
507 {
508   return vertices_.emplace(nil_point);
509 }
510 
new_vertex(const Point_d & p)511 Vertex_handle  new_vertex(const Point_d& p)
512 /*{\Mop adds a new vertex to |\Mvar| and returns it. The new vertex
513         has |p| as the associated point, but is has no associated
514         simplex nor index yet.}*/
515 {
516   return vertices_.emplace(p);
517 }
518 
associate_vertex_with_simplex(Simplex_handle s,int i,Vertex_handle v)519 void associate_vertex_with_simplex(Simplex_handle s, int i, Vertex_handle v)
520 /*{\Mop sets the $i$-th vertex of |s| to |v| and records this fact in
521 $v$. The latter occurs only if $v$ is non-nil.}*/
522 { s -> set_vertex(i,v);
523   if ( v != Vertex_handle() ) {
524     v -> set_simplex(s); v -> set_index(i);
525   }
526 }
527 
associate_point_with_vertex(Vertex_handle v,const Point_d & p)528 void associate_point_with_vertex(Vertex_handle v, const Point_d& p)
529 /*{\Mop sets the point associated with $v$ to $p$.}*/
530 { v -> set_point(p); }
531 
set_neighbor(Simplex_handle s,int i,Simplex_handle s1,int j)532 void set_neighbor(Simplex_handle s, int i, Simplex_handle s1, int j)
533 /*{\Mop sets the neihbor opposite to vertex $i$ of |s| to |s1| and
534         records vertex $j$ of |s1| as the vertex opposite to $i$.}*/
535 { s  -> set_neighbor(i,s1);
536   s1 -> set_neighbor(j,s);
537   s  -> set_opposite_vertex_index(i,j);
538   s1 -> set_opposite_vertex_index(j,i);
539 }
540 
541 void check_topology() const;
542 /*{\Mop Partially checks whether |\Mvar| is an abstract simplicial
543 complex. This function terminates without error if each vertex is a
544 vertex of the simplex of which it claims to be a vertex, if the
545 vertices of all simplices are pairwise distinct, if the neighbor
546 relationship is symmetric, and if neighboring simplices share exactly
547 |dcur| vertices.  It returns an error message if one of these
548 conditions is violated.  Note that it is not checked whether simplices
549 that share |dcur| vertices are neighbors in the data structure.}*/
550 
551 void check_topology_and_geometry() const;
552 /*{\Mop In addition to the above, this function checks whether all
553 vertices have an associated point different from
554 |Regular_complex_d::nil_point| and whether the points associated with the
555 vertices of any simplex are affinely independent. It returns an error
556 message otherwise.  Note that it is not checked whether the
557 intersection of any two simplices is a facet of both.}*/
558 
559 
560 typedef size_t Size_type;
561 
number_of_vertices()562 Size_type number_of_vertices() const  { return this->vertices_.size();}
number_of_simplices()563 Size_type number_of_simplices() const  { return this->simplices_.size();}
564 
565 void print_statistics(std::ostream& os = std::cout) const
566 {
567   os << "Regular_complex_d - statistic" << std::endl;
568   os << "number of vertices = " << number_of_vertices() << std::endl;
569   os << "number of simplices = " << number_of_simplices() << std::endl;
570 }
571 
572 /*{\Mtext \headerline{Lists and Iterators}
573 \setopdims{4.5cm}{3.5cm}}*/
574 
575 /*{\Mtext The following operation pairs return iterator ranges in
576 the style of STL.}*/
577 
vertices_begin()578 Vertex_iterator vertices_begin() { return vertices_.begin(); }
579 /*{\Mop the first vertex of |\Mvar|.}*/
vertices_end()580 Vertex_iterator vertices_end()   { return vertices_.end(); }
581 /*{\Mop the beyond vertex of |\Mvar|.}*/
simplices_begin()582 Simplex_iterator simplices_begin() { return simplices_.begin(); }
583 /*{\Mop the first simplex of |\Mvar|.}*/
simplices_end()584 Simplex_iterator simplices_end()   { return simplices_.end(); }
585 /*{\Mop the beyond simplex of |\Mvar|.}*/
586 
vertices_begin()587 Vertex_const_iterator vertices_begin() const { return vertices_.begin(); }
vertices_end()588 Vertex_const_iterator vertices_end()   const { return vertices_.end(); }
simplices_begin()589 Simplex_const_iterator simplices_begin() const { return simplices_.begin(); }
simplices_end()590 Simplex_const_iterator simplices_end()   const { return simplices_.end(); }
591 
all_simplices()592 std::list<Simplex_handle> all_simplices()
593 /*{\Mop returns the set of all maximal simplices in |\Mvar|.}*/
594 { std::list<Simplex_handle> res; Simplex_iterator it;
595   forall_rc_simplices(it,*this) res.push_back(it);
596   return res; }
597 
all_simplices()598 std::list<Simplex_const_handle> all_simplices() const
599 { std::list<Simplex_const_handle> res; Simplex_const_iterator it;
600   forall_rc_simplices(it,*this) res.push_back(it);
601   return res; }
602 
all_vertices()603 std::list<Vertex_handle> all_vertices()
604 /*{\Mop returns the set of all vertices in |\Mvar|.}*/
605 { std::list<Vertex_handle> res; Vertex_iterator it;
606   forall_rc_vertices(it,*this) res.push_back(it);
607   return res; }
608 
all_vertices()609 std::list<Vertex_const_handle> all_vertices() const
610 { std::list<Vertex_const_handle> res; Vertex_const_iterator it;
611   forall_rc_vertices(it,*this) res.push_back(it);
612   return res; }
613 
614 
615 
kernel()616 const R& kernel() const { return Kernel_; }
617 static const Point_d nil_point;
618 
619 }; // Regular_complex_d<R>
620 
621 // init static member:
622 template <class R>
623 const typename Regular_complex_d<R>::Point_d Regular_complex_d<R>::nil_point;
624 
625 
626 template <class R>
check_topology()627 void Regular_complex_d<R>::check_topology() const
628 {
629   Simplex_const_handle s,t;
630   Vertex_const_handle v;
631   int i,j,k;
632   if (dcur == -1) {
633     if (!vertices_.empty() || !simplices_.empty() )
634       CGAL_error_msg(      "check_topology: dcur is -1 but there are vertices or simplices");
635   }
636 
637   forall_rc_vertices(v,*this) {
638     if ( v != vertex(simplex(v),index(v)) )
639       CGAL_error_msg(      "check_topology: vertex-simplex relationship corrupted");
640   }
641 
642   forall_rc_simplices(s,*this) {
643     for(i = 0; i <= dcur; i++) {
644       for (j = i + 1; j <= dcur; j++) {
645         if (vertex(s,i) == vertex(s,j))
646           CGAL_error_msg(          "check_topology: a simplex with two equal vertices");
647       }
648     }
649   }
650 
651   forall_rc_simplices(s,*this) {
652     for(i = 0; i <= dcur; i++) {
653       if ((t = opposite_simplex(s,i)) != Simplex_const_handle()) {
654         int l = index_of_opposite_vertex(s,i);
655         if (s != opposite_simplex(t,l) ||
656             i != index_of_opposite_vertex(t,l))
657           CGAL_error_msg(          "check_topology: neighbor relation is not symmetric");
658 
659         for (j = 0; j <= dcur; j++) {
660           if (j != i) {
661             // j must also occur as a vertex of t
662             for (k = 0; k <= dcur &&
663                    ( vertex(s,j) != vertex(t,k) || k == l); k++) {}
664             if (k > dcur)
665               CGAL_error_msg(              "check_topology: too few shared vertices.");
666           }
667         }
668       }
669     }
670   }
671 }
672 
673 template <class R>
check_topology_and_geometry()674 void Regular_complex_d<R>::check_topology_and_geometry() const
675 {
676   check_topology();
677   Vertex_const_handle v;
678   forall_rc_vertices(v,*this) {
679     if ( v == Vertex_const_handle() ||
680          associated_point(v).identical(Regular_complex_d<R>::nil_point) )
681       CGAL_error_msg("check_topology_and_geometry: \
682       vertex with nil_point or no associated point.");
683   }
684 
685   typename R::Affinely_independent_d affinely_independent =
686     kernel().affinely_independent_d_object();
687   Simplex_const_handle s;
688   forall_rc_simplices(s,*this) {
689     std::vector<Point_d> A(dcur + 1);
690     for (int i = 0; i <= dcur; i++)
691       A[i] = associated_point(s,i);
692     if ( !affinely_independent(A.begin(),A.end()) )
693       CGAL_error_msg("check_topology_and_geometry: \
694       corners of some simplex are not affinely independent");
695   }
696 }
697 
698 
699 /*{\Mtext
700 \headerline{Iteration Statements}
701 
702 {\bf forall\_rc\_simplices}($s,C$)
703 $\{$ ``the simplices of $C$ are successively assigned to $s$'' $\}$
704 
705 {\bf forall\_rc\_vertices}($v,C$)
706 $\{$ ``the vertices of $C$ are successively assigned to $v$'' $\}$
707 
708 }*/
709 
710 /*{\Mimplementation Each simplex stores its vertices, the adjacent
711 simplices, and the opposite vertices in arrays. The space requirement
712 for a simplex is $3 * |dim| * 4$ Bytes for the contents of the arrays
713 plus the actual space for the points plus the constant space overhead
714 for the arrays (see the manual pages for arrays).
715 
716 The class |Regular_complex_d| needs constant space plus space for a list of
717 simplices (which is about 12 bytes per simplex). The total space
718 requirement is therefore about $12(|dim| + 2)$ bytes times the number
719 of simplices.  }*/
720 
721 
722 
723 } //namespace CGAL
724 #endif // CGAL_REGULAR_COMPLEX_D_H
725