File: Triangulation_segment_traverser_3.h

package info (click to toggle)
cgal 6.1.1-2
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 144,952 kB
  • sloc: cpp: 811,597; ansic: 208,576; sh: 493; python: 411; makefile: 286; javascript: 174
file content (1105 lines) | stat: -rw-r--r-- 39,408 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
// Copyright (c) 2012  INRIA Sophia-Antipolis (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
//
// $URL: https://github.com/CGAL/cgal/blob/v6.1.1/Triangulation_3/include/CGAL/Triangulation_segment_traverser_3.h $
// $Id: include/CGAL/Triangulation_segment_traverser_3.h 08b27d3db14 $
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
//
// Author(s): Thijs van Lankveld, Jane Tournois


#ifndef CGAL_TRIANGULATION_SEGMENT_TRAVERSER_3_H
#define CGAL_TRIANGULATION_SEGMENT_TRAVERSER_3_H

#include <CGAL/license/Triangulation_3.h>

#include <iostream>
#include <utility>

#include <CGAL/assertions.h>
#include <CGAL/Triangulation_utils_3.h>

#include <CGAL/Triangulation_data_structure_3.h>
#include <CGAL/Triangulation_cell_base_3.h>
#include <CGAL/Triangulation_vertex_base_3.h>
#include <CGAL/Triangulation_simplex_3.h>

#include <optional>

// If defined, type casting is done statically,
// reducing type-safety overhead.
#define CGAL_TST_ASSUME_CORRECT_TYPES

namespace CGAL {

template < class Tr, class Inc >
class Triangulation_segment_cell_iterator_3;

namespace internal {
template < class Tr >
struct Incrementer {
    typedef Incrementer<Tr>                                 Self;
    typedef Triangulation_segment_cell_iterator_3<Tr,Self>  SCI;    // describes the type of iterator expected by the incrementer.
    Incrementer() {}
    void increment( SCI& sci ) { sci.walk_to_next(); }
}; // struct Incrementer

} // namespace internal


//  provides an iterator over the cells intersected by a line segment.
/*
 *        The `Triangulation_segment_traverser_3` iterates over the cells
 *        of a `Triangulation_3` by following a straight line segment \f$ st \f$.
 *
 *        This class is closely related to `Triangulation_3::locate(...)`.
 *        However, unlike this `locate(...)` method, all the cells traversed
 *        by the `Triangulation_segment_traverser_3` intersect the interior of the line
 *        segment \f$ st \f$.
 *
 *        Traversal starts from a cell containing \f$ s \f$ and it ends in a cell containing
 *        \f$ t \f$.
 *        If \f$ st \f$ is coplanar with a facet or collinear with an edge, at most one of the
 *        incident cells is traversed.
 *        If \f$ st \f$ intersects an edge or vertex, at most two incident cells are traversed:
 *        the cells intersected by \f$ st \f$ strictly in their interior.
 *
 *        If \f$ s \f$ lies on the convex hull, traversal starts in an incident cell inside
 *        the convex hull. Similarly, if \f$ t \f$ lies on the convex hull, traversal ends in
 *        an adjacent cell inside the convex hull.
 *
 *        Both \f$ s \f$ and \f$ t \f$ may lie outside the convex hull of the triangulation,
 *        but they must lie within the affine hull of the triangulation. In either case, the
 *        finite facet of any infinite cells traversed must intersect \f$ st \f$.
 *
 *        The traverser may be applied to any triangulation of dimension > 0.
 *        However, for triangulations of dimension 1, the functionality is somewhat trivial.
 *
 *        The traverser becomes invalid whenever the triangulation is changed.
 *
 *        \tparam Tr_ is the triangulation type to traverse.
 *
 *        \cgalModels{ForwardIterator}
 *
 *        \sa `Triangulation_3`
 *        \sa `Forward_circulator_base`
 */
template < class Tr_, class Inc = internal::Incrementer<Tr_> >
class Triangulation_segment_cell_iterator_3
{
    typedef Tr_                                         Tr;
    typedef typename Tr::Triangulation_data_structure   Tds;
    typedef typename Tr::Geom_traits                    Gt;
    typedef Inc                                         Incrementer;

public:
// \name Types
// \{
    typedef Tr                                          Triangulation;          //< defines the triangulation type.
    typedef Triangulation_segment_cell_iterator_3<Tr,Inc>
                                                        Segment_cell_iterator;  //< defines the segment cell iterator type.

    typedef typename Tr::Point                          Point;                  //< defines the point type.
    typedef typename Tr::Segment                        Segment;                //< defines the line segment type.

    typedef typename Tr::Cell                           Cell;                   //< defines the type of a cell of the triangulation.
    typedef typename Tr::Edge                           Edge;                   //< defines the type of an edge of the triangulation.
    typedef typename Tr::Facet                          Facet;                  //< defines the type of a facet of the triangulation.

    typedef typename Tr::Vertex_handle                  Vertex_handle;          //< defines the type of a handle for a vertex in the triangulation.
    typedef typename Tr::Cell_handle                    Cell_handle;            //< defines the type of a handle for a cell in the triangulation.

    typedef typename Tr::Locate_type                    Locate_type;            //< defines the simplex type returned from location.

    struct Simplex                                                              //< defines the simplex type
    {
      Cell_handle cell = {};
      Locate_type lt = Locate_type::OUTSIDE_AFFINE_HULL;
      int li = -1;
      int lj = -1;
    };

    typedef Cell                                        value_type;             //< defines the value type the iterator refers to.
    typedef Cell&                                       reference;              //< defines the reference type of the iterator.
    typedef Cell*                                       pointer;                //< defines the pointer type of the iterator.
    typedef std::size_t                                 size_type;              //< defines the integral type that can hold the size of a sequence.
    typedef std::ptrdiff_t                              difference_type;        //< defines the signed integral type that can hold the distance between two iterators.
    typedef std::forward_iterator_tag                   iterator_category;      //< defines the iterator category.
// \}

    // describes the iterator type when applied to another type of triangulation or incrementer.
    template < class Tr2, class Inc2 >
    struct Rebind { typedef Triangulation_segment_cell_iterator_3<Tr2,Inc2>  Other; };

#if CGAL_DEBUG_TRIANGULATION_SEGMENT_TRAVERSER_3
    static auto display_vert(Vertex_handle v)
    {
      std::stringstream os;
      os.precision(17);
      if(v->time_stamp() == 0) {
        os << "inf";
      } else {
        os << '#' << v->time_stamp() << "=(" << v->point() << ")";
      }
      return os.str();
    };

    static auto display_lt(Locate_type lt) {
      std::stringstream os;
      switch(lt) {
        case Locate_type::VERTEX: os << " VERTEX"; break;
        case Locate_type::EDGE: os << " EDGE"; break;
        case Locate_type::FACET: os << " FACET"; break;
        case Locate_type::CELL: os << " CELL"; break;
        case Locate_type::OUTSIDE_CONVEX_HULL: os << " OUTSIDE_CONVEX_HULL"; break;
        case Locate_type::OUTSIDE_AFFINE_HULL: os << " OUTSIDE_AFFINE_HULL"; break;
      }
      return os.str();
    }

    static auto debug_simplex(Simplex s) {
      std::stringstream os;
      os.precision(17);
      const auto [c, lt, i, j] = s;
      if(c == Cell_handle{}) {
        os << "end()";
      } else {
        os << display_vert(c->vertex(0)) << " - " << display_vert(c->vertex(1)) << " - "
           << display_vert(c->vertex(2)) << " - " << display_vert(c->vertex(3));
        os << display_lt(lt) << " " << i << " " << j;
      }
      return os.str();
    }

    auto debug_iterator() const
    {
      std::stringstream os;
      os.precision(17);
      os << "  prev: " << debug_simplex(_prev) << "\n  cur: " << debug_simplex(_cur);
      return os.str();
    }
#endif // CGAL_DEBUG_TRIANGULATION_SEGMENT_TRAVERSER_3

private:
    typedef Segment_cell_iterator                       SCI;

    friend internal::Incrementer<Tr>;

protected:
// \internal \name Protected Attributes
// \{
    // \internal The triangulation to traverse.
    const Tr*       _tr;

// \}

    // The source and target points of the traversal.
    // These are also stored as vertices for cheaper equality computation.
    Point           _source;
    Point           _target;
    Vertex_handle   _s_vertex;
    Vertex_handle   _t_vertex;

    // The current cell with its entry point and the previous cell with its
    // exit point.
    // Note that the current cell will be Cell_handle() after incrementing past
    // the first cell containing the target.
    Simplex _cur, _prev;

public:
// \name Constructors
// \{
    //  constructs an iterator.
    /*  \param tr the triangulation to iterate though. This triangulation must have dimension > 0.
     *        \param s the source vertex. This vertex must be initialized and cannot be the infinite vertex.
     *        \param t the target vertex. This vertex must be initialized and cannot be the infinite vertex.
     *        It cannot equal `s`.
     */
        Triangulation_segment_cell_iterator_3( const Tr* tr, Vertex_handle s, Vertex_handle t );

    //  constructs an iterator.
    /*  \param tr the triangulation to iterate though. This triangulation must have dimension > 0.
     *        \param s the source vertex. This vertex must be initialized and cannot be the infinite vertex.
     *        \param t the target point. This point must be initialized and it cannot be be at the same location as `s`.
     *        If `tr` has dimension < 3, `t` must lie inside the affine hull of `tr`.
     */
        Triangulation_segment_cell_iterator_3( const Tr* tr, Vertex_handle s, const Point& t );

    //  constructs an iterator.
    /*  \param tr the triangulation to iterate though. This triangulation must have dimension > 0.
     *        \param s the source point. This point must be initialized and it cannot be be at the same location as `t`.
     *        \param t the target vertex. This vertex must be initialized and cannot be the infinite vertex.
     *        If `tr` has dimension < 3, `s` must lie inside the affine hull of `tr`.
     *        \param hint the starting point to search for `s`.
     */
        Triangulation_segment_cell_iterator_3( const Tr* tr, const Point& s, Vertex_handle t, Cell_handle hint = Cell_handle() );

    //  constructs an iterator.
    /*  \param tr the triangulation to iterate though. This triangulation must have dimension > 0.
     *        \param s the source point. This point must be initialized. If `tr` has dimension < 3, `s` must lie inside
     *        the affine hull of `tr`.
     *        \param t the target point. This point must be initialized and it cannot be be at the same location as `s`.
     *        If `tr` has dimension < 3, `t` must lie inside the affine hull of `tr`.
     *        \param hint the starting point to search for `s`.
     */
    Triangulation_segment_cell_iterator_3( const Tr* tr, const Point& s, const Point& t, Cell_handle hint = Cell_handle() );

    //  constructs an iterator.
    /*  \param tr the triangulation to iterate though. This triangulation must have dimension > 0.
     *        \param S the segment to be traversed. If `tr` has dimension < 3, `S` must lie inside
     *        the affine hull of `tr`. `S` must not be degenerate, i.e. its source and target must not be equal.
     *        \param hint the starting point to search for `S`.
     */
    Triangulation_segment_cell_iterator_3( const Tr* tr, const Segment& S, Cell_handle hint = Cell_handle() );
// \}

    // private constructor that does not initialize the source and target.
    // used for the end()
    Triangulation_segment_cell_iterator_3(const Tr* tr);

#ifndef CGAL_TST_ASSUME_CORRECT_TYPES
    // The virtual destructor is mainly defined to indicate to the casting
    // operators that this is a dynamic type.
    virtual
#endif
    ~Triangulation_segment_cell_iterator_3() {}


public:
// \name Accessors
// \{

    const Tr* triangulation() const     { return _tr; }

    //  gives the source point of the segment followed.
    /*  \return the source point.
     */
    const Point&    source() const      { return _source; }

    //  gives the target point of the segment followed.
    /*  \return the target point.
         */
    const Point&    target() const      { return _target; }

    Vertex_handle target_vertex() const { return _t_vertex; }

    //  gives a handle to the current cell.
    /*  By invariance, this cell is intersected by the segment
         *        between `source()` and `target()`.
         *        \return a handle to the current cell.
         *        \sa `cell()`.
         */
    Cell_handle     handle() const
    {
      return _cur.cell;
    }

    //  gives the previous cell.
    /*  This cell is uninitialized until the iterator leaves the initial
         *        cell.
         *        By invariance, once initialized, this cell must be intersected by the segment
         *        between `source()` and `target()`.
         *        \return the previous cell.
         *        \sa `handle()`.
         */
    Cell_handle     previous() const
    {
      return prev_cell();
    }

    //  provides a dereference operator.
    /*         \return a pointer to the current cell.
         */
    Cell*           operator->()
    {
      return &*(_cur.cell);
    }

    //  provides an indirection operator.
    /*  \return the current cell.
         */
    Cell&           operator*()
    {
      return *(_cur.cell);
    }

    //  provides a conversion operator.
    /*         \return a handle to the current cell.
         */
    operator const Cell_handle&() const
    {
      return _cur.cell;
    }

    //  provides a conversion operator.
    /*         \return the simplex through which the current cell was entered.
         */
    operator Simplex() const { return _cur; }

    //  checks whether the iterator has reached the final cell, which contains the `target()`.
    /*  If the `target()` lies on a facet, edge, or vertex, the final cell is the cell containing
         *        the interior of the segment between `source()` and `target()`.
         *        \return true iff the current cell contains the `target()`.
         */
    bool            has_next() const
    {
      return this->cell() != Cell_handle();
    }

    //  gives the simplex through which the current cell was entered.
    /*         For the first cell, containing the `source()` \f$ s \f$,
     *  this indicates the location of \f$ s \f$ in this cell.
         */
    void            entry( Locate_type& lt, int& li, int& lj ) const
    {
      lt = this->lt(); li = this->li(); lj = this->lj();
    }
    std::tuple<Locate_type, int, int> entry() const
    {
      return { lt(), li(), lj() };
    }
    //  gives the simplex through which the previous cell was exited.
    /*         \pre the current cell is not the initial cell.
         */
    void            exit( Locate_type& lt, int& li, int& lj ) const
    {
      lt = prev_lt(); li = prev_li(); lj = prev_lj();
    }
    std::tuple<Locate_type, int, int> exit() const
    {
      return { prev_lt(), prev_li(), prev_lj() };
    }

    //  gives the past-the-end iterator associated with this iterator.
    SCI             end() const;
// \}

public:
// \name Mutators
// \{
    //  provides the increment postfix operator.
    /*         After incrementing the iterator, the current cell intersects the segment
     *        between `source()` and `target()` closer to the `target()` than the previous cell.
     *        \sa `operator++(int)`.
     *  \pre The current cell does not contain the `target()`.
     */
    SCI&            operator++();

    //  provides the increment prefix operator.
    /*         After incrementing the iterator, the current cell intersects the segment
     *        between `source()` and `target()` closer to the `target()` than the previous cell.
     *        than the previous cell.
     *        \sa `operator++()`.
     *  \pre The current cell does not contain the `target()`.
     */
    SCI             operator++( int );

    //  iterates to the final cell, which contains the `target()`.
    /*         \return the final cell.
     */
    Cell_handle     complete();
// \}

public:
// \name Comparison
// \{
    //  compares this iterator with `sci`.
    /*  \param sci the other iterator.
     *        \return true iff the other iterator iterates the same triangulation along the same line segment
     *        and has the same current cell.
     *        \sa `operator!=( const SCI& t )`.
     */
    bool            operator==( const SCI& sci ) const;

    //  compares this iterator with `sci`.
    /*  \param sci the other iterator.
     *        \return `false` iff the other iterator iterates the same triangulation along the same line segment
     *        and has the same current cell.
     *        \sa `operator==( const SCI& t ) const`.
     */
    bool            operator!=( const SCI& sci ) const;

    //  compares the current cell with `ch`.
    /*  \param ch a handle to the other cell.
     *        \return true iff the current cell is the same as the one pointed to by `ch`.
     *        \sa `operator!=( const Cell_handle& ch ) const`.
     *        \sa `operator==( typename TriangulationTraits_3::Cell_handle ch, Triangulation_segment_cell_iterator_3<TriangulationTraits_3> t )`.
     */
    bool            operator==( const Cell_handle& ch ) const
    {
      return ch == _cur.cell;
    }

    //  compares the current cell with `ch`.
    /*  \param ch a handle to the other cell.
     *        \return `false` iff the current cell is the same as the one pointed to by `ch`.
     *        \sa `operator==( const Cell_handle& ch )`.
     *        \sa `operator!=( typename TriangulationTraits_3::Cell_handle ch, Triangulation_segment_cell_iterator_3<TriangulationTraits_3> t )`.
     */
    bool            operator!=( const Cell_handle& ch ) const
    {
      return ch != _cur.cell;
    }
// \}

        bool            operator==( Nullptr_t CGAL_assertion_code(n) ) const;
        bool            operator!=( Nullptr_t n ) const;

protected:
// \internal \name Protected Member Functions
// \{
    //  walks to the next cell.
    /*  \sa `complete()`.
         */
    void            walk_to_next();

    //  increments the iterator.
    /*  This method may perform more actions based on the superclass.
     *  \sa `complete()`.
         */
    void            increment() {
        typedef typename Incrementer::SCI    Expected;
#ifdef CGAL_TST_ASSUME_CORRECT_TYPES
        Expected& sci = static_cast<Expected&>( *this );
#else // CGAL_TST_ASSUME_CORRECT_TYPES
        Expected& sci = dynamic_cast<Expected&>( *this );
#endif // CGAL_TST_ASSUME_CORRECT_TYPES
        Incrementer().increment( sci );
    }
// \}

private:
    // at the end of the constructors, entry() is a vertex, edge or facet,
    // we need to circulate/iterate over its incident cells to
    // make sure that the current cell intersects the input query
    void jump_to_intersecting_cell();

    //  walk_to_next(), if the triangulation is 3D.
    std::pair<Simplex, Simplex> walk_to_next_3(const Simplex& prev,
                                               const Simplex& cur) const;
    void            walk_to_next_3_inf( int inf );

    //  walk_to_next(), if the triangulation is 2D.
    void            walk_to_next_2();
    void            walk_to_next_2_inf( int inf );

private:
    inline int      edgeIndex( int i, int j ) const {
        CGAL_precondition( i>=0 && i<=3 );
        CGAL_precondition( j>=0 && j<=3 );
        CGAL_precondition( i != j );
        return ( i==0 || j==0 ) ? i+j-1 : i+j;
    }

    bool have_same_entry(const Simplex& s1, const Simplex& s2) const;

    // Compute the orientation of a point compared to the oriented plane supporting a half-facet.
    CGAL::Orientation orientation(const Facet& f, const Point& p) const;

    bool coplanar(const Facet &f, const Point &p) const;

    // Gives the edge incident to the same cell that is not incident to any of the input vertices.
    Edge opposite_edge(Cell_handle c, int li, int lj) const;
    Edge opposite_edge(const Edge& e) const;

protected:
    // ref-accessors to the simplex, for use in internal code
    // access _cur
    Cell_handle& cell()             { return _cur.cell; }
    Cell_handle const& cell() const { return _cur.cell; }

    Locate_type& lt()             { return _cur.lt; }
    Locate_type const& lt() const { return _cur.lt; }

    int& li()             { return _cur.li; }
    int const& li() const { return _cur.li; }

    int& lj()             { return _cur.lj; }
    int const& lj() const { return _cur.lj; }

    // access _prev
    Cell_handle& prev_cell()             { return _prev.cell; }
    Cell_handle const& prev_cell() const { return _prev.cell; }

    Locate_type& prev_lt()             { return _prev.lt; }
    Locate_type const& prev_lt() const { return _prev.lt; }

    int& prev_li()             { return _prev.li; }
    int const& prev_li() const { return _prev.li; }

    int& prev_lj()             { return _prev.lj; }
    int const& prev_lj() const { return _prev.lj; }

}; // class Triangulation_segment_cell_iterator_3

//  compares a handle to a cell to a traverser.
/*  \param ch the handle to a cell.
 *        \param t the traverser.
 *        \return true iff the cell currently traversed by `t` is the same as the one pointed to by `ch`.
 *        \sa `operator!=( typename TriangulationTraits_3::Cell_handle ch, Triangulation_segment_cell_iterator_3<TriangulationTraits_3> t )`.
 *        \sa `Triangulation_segment_cell_iterator_3::operator==( const Cell_handle& ch )`.
 */
template < class Tr, class Inc >
inline bool operator==( typename Tr::Cell_handle ch, Triangulation_segment_cell_iterator_3<Tr,Inc> tci ) { return tci == ch; }

//  compares a handle to a cell to a traverser.
/*  \param ch the handle to a cell.
 *        \param t the traverser.
 *        \return `false` iff the cell currently traversed by `t` is the same as the one pointed to by `ch`.
 *        \sa `operator==( typename TriangulationTraits_3::Cell_handle ch, Triangulation_segment_cell_iterator_3<TriangulationTraits_3> t )`.
 *        \sa `Triangulation_segment_cell_iterator_3::operator!=( const Cell_handle& ch )`.
 */
template < class Tr, class Inc >
inline bool operator!=( typename Tr::Cell_handle ch, Triangulation_segment_cell_iterator_3<Tr,Inc> tci ) { return tci != ch; }



/********************************************************************/
/********************************************************************/
/********************************************************************/
template < class Tr_, class Inc = internal::Incrementer<Tr_> >
class Triangulation_segment_simplex_iterator_3
{
  typedef Tr_                                         Tr;
  typedef typename Tr::Triangulation_data_structure   Tds;
  typedef typename Tr::Geom_traits                    Gt;
  typedef Inc                                         Incrementer;

private:
  typedef Triangulation_segment_simplex_iterator_3<Tr_, Inc> Simplex_iterator;
  typedef Triangulation_segment_cell_iterator_3<Tr_, Inc>    SCI;

private:
  typedef typename SCI::Point    Point;
  typedef typename SCI::Segment  Segment;

public:
  // \{
  typedef typename SCI::Vertex_handle Vertex_handle;//< defines the type of a handle for a vertex in the triangulation
  typedef typename SCI::Cell_handle   Cell_handle;  //< defines the type of a handle for a cell in the triangulation.
  typedef typename SCI::Cell          Cell;  //< defines the type of a handle for a cell in the triangulation.
  typedef typename SCI::Triangulation::Edge  Edge;  //< defines the type of an edge in the triangulation.
  typedef typename SCI::Triangulation::Facet Facet; //< defines the type of a facet in the triangulation.
  typedef typename SCI::Locate_type   Locate_type;  //< defines the simplex type returned from location.

  typedef CGAL::Triangulation_simplex_3<Tds>  Simplex_3;

  typedef Simplex_3   value_type;       //< defines the value type  the iterator refers to.
  typedef const Simplex_3&  reference;  //< defines the reference type of the iterator.
  typedef const Simplex_3*  pointer;    //< defines the pointer type of the iterator.
  typedef std::size_t    size_type;        //< defines the integral type that can hold the size of a sequence.
  typedef std::ptrdiff_t difference_type;  //< defines the signed integral type that can hold the distance between two iterators.
  typedef std::forward_iterator_tag iterator_category;      //< defines the iterator category.
  // \}

private:
  SCI _cell_iterator;
  Simplex_3 _curr_simplex;

public:
  Triangulation_segment_simplex_iterator_3(const Tr* tr
    , Vertex_handle s, Vertex_handle t)
    : _cell_iterator(tr, s, t)
  { set_curr_simplex_to_entry(); }
  Triangulation_segment_simplex_iterator_3(const Tr* tr
    , Vertex_handle s, const Point& t)
    : _cell_iterator(tr, s, t)
  { set_curr_simplex_to_entry(); }
  Triangulation_segment_simplex_iterator_3(const Tr* tr
    , const Point& s, Vertex_handle t, Cell_handle hint = Cell_handle())
    : _cell_iterator(tr, s, t, hint)
  { set_curr_simplex_to_entry(); }
  Triangulation_segment_simplex_iterator_3(const Tr* tr
    , const Point& s, const Point& t, Cell_handle hint = Cell_handle())
    : _cell_iterator(tr, s, t, hint)
  { set_curr_simplex_to_entry(); }
  Triangulation_segment_simplex_iterator_3(const Tr* tr
    , const Segment& seg, Cell_handle hint = Cell_handle())
    : _cell_iterator(tr, seg, hint)
  { set_curr_simplex_to_entry(); }
  Triangulation_segment_simplex_iterator_3(const Tr* tr)
    : _cell_iterator(tr)
    , _curr_simplex()
  {}

  bool operator==(const Simplex_iterator& sit) const
  {
    return sit._cell_iterator == _cell_iterator
        && sit._curr_simplex == _curr_simplex;
  }
  bool operator!=(const Simplex_iterator& sit) const
  {
    return sit._cell_iterator != _cell_iterator
        || sit._curr_simplex != _curr_simplex;
  }

  const Point&    source() const      { return _cell_iterator.source(); }
  const Point&    target() const      { return _cell_iterator.target(); }

  const Tr& triangulation() const     { return *_cell_iterator.triangulation(); }

private:
  Triangulation_segment_simplex_iterator_3
    (const SCI& sci)
    : _cell_iterator(sci)
    , _curr_simplex()
  {}

private:
  void set_curr_simplex_to_entry()
  {
#if CGAL_DEBUG_TRIANGULATION_SEGMENT_TRAVERSER_3
    std::cerr << "cell iterator is:\n" << _cell_iterator.debug_iterator() << std::endl;
#endif // #if CGAL_DEBUG_TRIANGULATION_SEGMENT_TRAVERSER_3

    Locate_type lt;
    int li, lj;
    Cell_handle cell = Cell_handle(_cell_iterator);

    //check what is the entry type of _cell_iterator
    if (cell == Cell_handle())
    {
      //where did the segment get out from previous cell
      cell = _cell_iterator.previous();
      _cell_iterator.exit(lt, li, lj);
    }
    else
    {
      _cell_iterator.entry(lt, li, lj);
    }

    switch (lt)
    {
    case Locate_type::VERTEX:
      _curr_simplex = cell->vertex(li);
      break;
    case Locate_type::EDGE:
      _curr_simplex = Edge(cell, li, lj);
      break;
    case Locate_type::FACET:
      _curr_simplex = Facet(cell, li);
      break;
      //the 3 cases below correspond to the case when _cell_iterator
      //is in its initial position: _cur is locate(source)
    case Locate_type::CELL:
    case Locate_type::OUTSIDE_CONVEX_HULL:
    case Locate_type::OUTSIDE_AFFINE_HULL:
      if (Cell_handle(_cell_iterator) == Cell_handle())
        _curr_simplex = Simplex_3();
      else
        _curr_simplex = cell;
      break;
    default:
      CGAL_unreachable();
    };
  }

public:
  Simplex_iterator end() const
  {
    Simplex_iterator sit(_cell_iterator.end());
    return sit;
  }

  //  provides the increment postfix operator.
  Simplex_iterator& operator++()
  {
    auto increment_cell_iterator = [&]() {
      ++_cell_iterator;
#if CGAL_DEBUG_TRIANGULATION_SEGMENT_TRAVERSER_3
      std::cerr << "increment cell iterator to:\n" << _cell_iterator.debug_iterator() << '\n';
#endif
    };
    CGAL_assertion(_curr_simplex.incident_cell() != Cell_handle());

    if(!cell_iterator_is_ahead()) {
      increment_cell_iterator(); // cell_iterator needs to be ahead
    }

    Cell_handle ch_next = Cell_handle(_cell_iterator);
    Cell_handle ch_prev = _cell_iterator.previous();
    Locate_type lt_prev;
    int li_prev, lj_prev;
    _cell_iterator.exit(lt_prev, li_prev, lj_prev);

    if(_curr_simplex.dimension() == 3) {
      set_curr_simplex_to_entry();
      return *this;
    }
    if(lt_prev == Locate_type::CELL ||
       lt_prev == Locate_type::OUTSIDE_CONVEX_HULL ||
       lt_prev == Locate_type::OUTSIDE_AFFINE_HULL)
    {
      CGAL_assertion(ch_next == Cell_handle());
      _curr_simplex = ch_prev;
      return *this;
    }

    switch(_curr_simplex.dimension()) {
    case 2: { /*Facet*/
      CGAL_assertion((ch_next == Cell_handle()) == (_cell_iterator == _cell_iterator.end()));

      switch(lt_prev) {
      case Locate_type::VERTEX: { // facet-cell?-vertex-outside
        Vertex_handle v_prev{ch_prev->vertex(li_prev)};
        if(facet_has_vertex(get_facet(), v_prev))
          _curr_simplex = v_prev;
        else
          _curr_simplex = ch_prev;
      } break;
      case Locate_type::EDGE: { // facet-cell?-edge-outside
        Edge edge_prev{ch_prev, li_prev, lj_prev};
        if(facet_has_edge(get_facet(), edge_prev))
          _curr_simplex = edge_prev;
        else
          _curr_simplex = ch_prev;
      } break;
      case Locate_type::FACET: { // facet-cell-facet-outside
        Facet f_prev{ch_prev, li_prev};
        if(is_same_facet(f_prev, get_facet())) {
          if(ch_next == Cell_handle())
            _curr_simplex = Simplex_3();
          else
            _curr_simplex = ch_next;
        } else
          _curr_simplex = ch_prev;
      } break;
      default:
        CGAL_unreachable();
      }
    } break;
    case 1: {/*Edge*/
      switch(lt_prev) {
      case Locate_type::VERTEX: { //edge-vertex-outside
        Vertex_handle v_prev{ch_prev->vertex(li_prev)};
        if(edge_has_vertex(get_edge(), v_prev))
          _curr_simplex = v_prev;
        else
          _curr_simplex = shared_facet(get_edge(), v_prev);
      } break;
      case Locate_type::EDGE: { //edge-outside or edge-cell-edge-outside
        const Edge e_prev(ch_prev, li_prev, lj_prev);
        if(is_same_edge(get_edge(), e_prev)) {
          if(ch_next == Cell_handle()) {
            _curr_simplex = Simplex_3();
          } else {
            _curr_simplex = ch_next;
          }
        } else {
          auto facet_opt = shared_facet(get_edge(), e_prev);
          if(facet_opt.has_value()) {
            _curr_simplex = facet_opt.value();
          }
          else {
            _curr_simplex = shared_cell(get_edge(), e_prev);
          }
        }
      } break;
      case Locate_type::FACET: {
        Facet f_prev{ch_prev, li_prev};
        if(facet_has_edge(f_prev, get_edge()))
          _curr_simplex = f_prev; //edge-facet-outside
        else
          _curr_simplex = ch_prev; //query goes through the cell
      } break;
      default:
        CGAL_unreachable();
      }
    } break;
    case 0 :/*Vertex_handle*/
    {
      switch(lt_prev) {
      case Locate_type::VERTEX: {
        if(ch_prev->vertex(li_prev) != get_vertex()) // avoid infinite loop edge-vertex-same edge-...
          _curr_simplex = Edge(ch_prev, li_prev, ch_prev->index(get_vertex()));
        else {
          if(ch_next == Cell_handle()) {
            _curr_simplex = Simplex_3();
          } else {
            _curr_simplex = ch_next;
          }
        }
      } break;
      case Locate_type::EDGE: {
        const Edge e_prev(ch_prev, li_prev, lj_prev);
        if(edge_has_vertex(e_prev, get_vertex()))
          _curr_simplex = e_prev;
        else
          _curr_simplex = shared_facet(Edge(ch_prev, li_prev, lj_prev), get_vertex());
      } break;
      case Locate_type::FACET: {
        if(ch_prev->vertex(li_prev) != get_vertex()) // vertex-facet-outside
          _curr_simplex = Facet(ch_prev, li_prev);
        else // vertex-cell-facet-outside
          _curr_simplex = ch_prev;
      } break;
      default:
        CGAL_unreachable();
      }
    } break;
    default:
      CGAL_unreachable();
    };
    return *this;
  }

  //  provides the increment prefix operator.
  Simplex_iterator operator++(int)
  {
    Simplex_iterator tmp(*this);
    ++(*this);
    return tmp;
  }

  //  provides a dereference operator.
  /*         \return a pointer to the current cell.
  */
  const Simplex_3*   operator->()        { return &_curr_simplex; }

  //  provides an indirection operator.
  /*  \return the current cell.
  */
  const Simplex_3&   operator*()         { return _curr_simplex; }

  //  provides a conversion operator.
  /*  \return the current simplex
  */
  operator const Simplex_3() const { return _curr_simplex; }

  bool is_vertex() const { return _curr_simplex.dimension() == 0; }
  bool is_edge()   const { return _curr_simplex.dimension() == 1; }
  bool is_facet()  const { return _curr_simplex.dimension() == 2; }
  bool is_cell()   const { return _curr_simplex.dimension() == 3; }

  const Cell cell() const
  {
    return _cell_iterator.cell();
  }

  const Simplex_3& get_simplex() const { return _curr_simplex; }
  Vertex_handle get_vertex() const
  {
    CGAL_assertion(is_vertex());
    return Vertex_handle(_curr_simplex);
  }
  Edge get_edge() const
  {
    CGAL_assertion(is_edge());
    return Edge(_curr_simplex);
  }
  Facet get_facet() const
  {
    CGAL_assertion(is_facet());
    return Facet(_curr_simplex);
  }
  Cell_handle get_cell() const
  {
    CGAL_assertion(is_cell());
    return Cell_handle(_curr_simplex);
  }

public:
  //returns true in any of the degenerate cases,
  //i.e. when _curr_simplex has the following values successively
  // edge / facet / edge
  // edge / facet / vertex
  // vertex / facet / edge
  // vertex / edge / vertex
  // TODO : rename this function
  bool is_collinear() const
  {
    int curr_dim = _curr_simplex.dimension();
    //this concerns only edges and facets
    if (curr_dim == 1 || curr_dim == 2)
      return cell_iterator_is_ahead();
      //the degeneracy has been detected by moving cell_iterator forward
    else
      return false;
  }

  int simplex_dimension() const
  {
    return _curr_simplex.dimension();
  }

private:
  bool cell_iterator_is_ahead() const
  {
    Cell_handle ch = Cell_handle(_cell_iterator);
    if(ch == Cell_handle())
      return true;

    switch (_curr_simplex.dimension())
    {
    case 0 ://vertex
      return !ch->has_vertex(get_vertex());
    case 1 ://edge
      return !cell_has_edge(ch, get_edge());
    case 2 ://facet
      return !cell_has_facet(ch, get_facet());
    case 3 ://cell
      return ch != get_cell();
    default:
      CGAL_unreachable();
    }
    //should not be reached
    CGAL_unreachable();
    return false;
  }

  bool cell_has_edge(const Cell_handle ch, const Edge& e) const
  {
    Vertex_handle v1 = e.first->vertex(e.second);
    Vertex_handle v2 = e.first->vertex(e.third);
    return ch->has_vertex(v1) && ch->has_vertex(v2);
  }
  bool cell_has_facet(const Cell_handle c, const Facet& f) const
  {
    return f.first == c
        || f.first->neighbor(f.second) == c;
  }

  bool facet_has_edge(const Facet& f, const Edge& e) const
  {
    Vertex_handle v1 = e.first->vertex(e.second);
    Vertex_handle v2 = e.first->vertex(e.third);
    Cell_handle c = f.first;
    const int fi = f.second;

    unsigned int count = 0;
    for (int i = 1; i < 4; ++i)
    {
      Vertex_handle vi = c->vertex((fi + i) % 4);
      if (vi == v1 || vi == v2)
        ++count;
      if (count == 2)
        return true;
    }
    return false;
  }

  bool facet_has_vertex(const Facet& f, const Vertex_handle v) const
  {
    return triangulation().tds().has_vertex(f, v);
  }

  bool edge_has_vertex(const Edge& e, const Vertex_handle v) const
  {
    return e.first->vertex(e.second) == v
        || e.first->vertex(e.third) == v;
  }

  bool is_same_edge(const Edge& e1, const Edge& e2) const
  {
    return edge_has_vertex(e1, e2.first->vertex(e2.second))
        && edge_has_vertex(e1, e2.first->vertex(e2.third));
  }

  bool is_same_facet(const Facet& f1, const Facet& f2) const
  {
    return f1 == f2 || triangulation().mirror_facet(f1) == f2;
  }

  std::optional<Vertex_handle> shared_vertex(const Edge& e1, const Edge& e2) const
  {
    Vertex_handle v1a = e1.first->vertex(e1.second);
    Vertex_handle v1b = e1.first->vertex(e1.third);
    Vertex_handle v2a = e2.first->vertex(e2.second);
    Vertex_handle v2b = e2.first->vertex(e2.third);

    if (v1a == v2a || v1a == v2b)
      return v1a;
    else if (v1b == v2a || v1b == v2b)
      return v1b;
    else
      return {};
  }

  std::optional<Facet> shared_facet(const Edge& e1, const Edge& e2) const
  {
    Vertex_handle v2a = e2.first->vertex(e2.second);
    Vertex_handle v2b = e2.first->vertex(e2.third);

    auto sv_opt = shared_vertex(e1, e2);
    if(!sv_opt.has_value())
      return {};
    Vertex_handle sv = sv_opt.value();
    Vertex_handle nsv2 = (sv == v2a) ? v2b : v2a;

    typename Tr::Facet_circulator circ
      = triangulation().incident_facets(e1);
    typename Tr::Facet_circulator end = circ;
    do
    {
      Facet f = *circ;
      for (int i = 1; i < 4; ++i)
      {
        if (nsv2 == f.first->vertex((f.second + i) % 4))
          return f;
      }
    } while (++circ != end);

    return {};
  }

  Facet shared_facet(const Edge& e, const Vertex_handle v) const
  {
    typename Tr::Facet_circulator circ
      = triangulation().incident_facets(e);
    typename Tr::Facet_circulator end = circ;
    do
    {
      Facet f = *circ;
      if (facet_has_vertex(f, v))
        return f;
    } while (++circ != end);

    std::cerr << "There is no facet shared by e and v" << std::endl;
    CGAL_unreachable();
    return Facet(Cell_handle(), 0);
  }

  Cell_handle shared_cell(const Edge& e, const Vertex_handle v) const
  {
    typename Tr::Cell_circulator circ
      = triangulation().incident_cells(e);
    typename Tr::Cell_circulator end = circ;
    do
    {
      Cell_handle c = circ;
      if (c->has_vertex(v))
        return c;
    } while (++circ != end);

    std::cerr << "There is no cell shared by e and v" << std::endl;
    CGAL_unreachable();
    return Cell_handle();
  }

  Cell_handle shared_cell(const Facet& f, const Vertex_handle v) const
  {
    Cell_handle c = f.first;
    if (c->has_vertex(v))
      return c;
    else
    {
      c = f.first->neighbor(f.second);
      CGAL_assertion(c->has_vertex(v));
      return c;
    }
  }

  Cell_handle shared_cell(const Edge e1, const Edge e2) const {
    Facet facet = shared_facet(e1, e2.first->vertex(e2.second));
    return shared_cell(facet, e2.first->vertex(e2.third));
  }

};//class Triangulation_segment_simplex_iterator_3

} // namespace CGAL

#include <CGAL/Triangulation_3/internal/Triangulation_segment_traverser_3_impl.h>

#endif // CGAL_TRIANGULATION_SEGMENT_TRAVERSER_3_H