File: clip.h

package info (click to toggle)
cgal 6.1.1-1
  • links: PTS, VCS
  • area: main
  • in suites: forky
  • size: 144,952 kB
  • sloc: cpp: 811,597; ansic: 208,576; sh: 493; python: 411; makefile: 286; javascript: 174
file content (1502 lines) | stat: -rw-r--r-- 61,512 bytes parent folder | download | duplicates (2)
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
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
// Copyright (c) 2016-2025 GeometryFactory (France).
// All rights reserved.
//
// This file is part of CGAL (www.cgal.org).
//
// $URL: https://github.com/CGAL/cgal/blob/v6.1.1/Polygon_mesh_processing/include/CGAL/Polygon_mesh_processing/clip.h $
// $Id: include/CGAL/Polygon_mesh_processing/clip.h 08b27d3db14 $
// SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial
//
//
// Author(s)     : Sebastien Loriot

#ifndef CGAL_POLYGON_MESH_PROCESSING_CLIP_H
#define CGAL_POLYGON_MESH_PROCESSING_CLIP_H

#include <CGAL/license/Polygon_mesh_processing/corefinement.h>

#include <CGAL/Polygon_mesh_processing/corefinement.h>
#include <CGAL/Polygon_mesh_processing/connected_components.h>
#include <CGAL/Polygon_mesh_processing/bbox.h>
#include <CGAL/Polygon_mesh_processing/triangulate_faces.h>
#include <CGAL/Polygon_mesh_processing/orientation.h>
#include <CGAL/Polygon_mesh_processing/triangulate_hole.h>
#include <CGAL/Polygon_mesh_processing/border.h>
#include <CGAL/Polygon_mesh_processing/repair.h>
#include <CGAL/Polygon_mesh_processing/internal/Corefinement/Generic_clip_output_builder.h>
#include <CGAL/iterator.h>

#include <CGAL/AABB_triangle_primitive_3.h>

#include <CGAL/boost/graph/properties.h>
#include <CGAL/boost/graph/Face_filtered_graph.h>

#include <boost/property_map/property_map.hpp>

#include <array>
#include <algorithm>
#include <map>
#include <set>
#include <type_traits>
#include <unordered_map>
#include <utility>
#include <vector>
#include <bitset>

#include <CGAL/Polygon_mesh_processing/refine_with_plane.h>

#ifndef CGAL_PLANE_CLIP_DO_NOT_USE_TRIANGULATION
#include <CGAL/Triangulation_vertex_base_with_info_2.h>
#include <CGAL/Triangulation_face_base_with_info_2.h>
#include <CGAL/Constrained_Delaunay_triangulation_2.h>
#include <CGAL/Projection_traits_3.h>
#include <CGAL/mark_domain_in_triangulation.h>
#include <CGAL/Delaunay_mesh_face_base_2.h>
#include <boost/iterator/transform_iterator.hpp>
#endif

namespace CGAL {
namespace Polygon_mesh_processing {
namespace internal {

#ifndef CGAL_PLANE_CLIP_DO_NOT_USE_TRIANGULATION
template <class Traits,
          class TriangleMesh,
          class VertexPointMap,
          class Visitor>
void close_and_triangulate(TriangleMesh& tm, VertexPointMap vpm, typename Traits::Vector_3 plane_normal, Visitor& visitor)
{
  using vertex_descriptor = typename boost::graph_traits<TriangleMesh>::vertex_descriptor;
  using halfedge_descriptor = typename boost::graph_traits<TriangleMesh>::halfedge_descriptor;
  using face_descriptor = typename boost::graph_traits<TriangleMesh>::face_descriptor;
  using P_traits = Projection_traits_3<Traits>;
  using Vb = Triangulation_vertex_base_with_info_2<vertex_descriptor, P_traits>;
  using Fb = CGAL::Delaunay_mesh_face_base_2<P_traits>;
  using TDS = Triangulation_data_structure_2<Vb, Fb>;
  using Itag = No_constraint_intersection_tag;
  using CDT = Constrained_Delaunay_triangulation_2<P_traits, TDS, Itag>;

  P_traits p_traits(plane_normal);

  std::vector<std::pair<typename Traits::Point_3, vertex_descriptor>> points;
  std::vector<std::pair<vertex_descriptor, vertex_descriptor>> csts;

  halfedge_descriptor href = boost::graph_traits<TriangleMesh>::null_halfedge ();
  bool first=true;

  for (halfedge_descriptor h : halfedges(tm))
  {
    if (is_border(h, tm))
    {
      if (first)
      {
        href=h;
        first=false;
      }
      vertex_descriptor src = source(h, tm), tgt = target(h, tm);
      points.emplace_back(get(vpm, tgt), tgt);
      csts.emplace_back(src, tgt);
    }
  }

  if (points.empty()) return;

  CDT cdt(p_traits);
  cdt.insert(points.begin(), points.end());

  // TODO: in case of degenerate edges, we don't triangulate but it's kind of an issue on our side
  if (cdt.number_of_vertices()!=points.size()) return;


  typedef CGAL::dynamic_vertex_property_t<typename CDT::Vertex_handle> V_tag;
  typename boost::property_map<TriangleMesh, V_tag>::type v2v = get(V_tag(), tm);
  for (auto vh : cdt.finite_vertex_handles())
    put(v2v, vh->info(), vh);

  for (auto vv : csts)
    cdt.insert_constraint(get(v2v, vv.first), get(v2v, vv.second));

  mark_domain_in_triangulation(cdt);

  typename CDT::Edge edge_ref;
  CGAL_assertion_code(bool OK =)
  cdt.is_edge(get(v2v, csts[0].first), get(v2v, csts[0].second), edge_ref.first, edge_ref.second);
  CGAL_assertion(OK);

  if (!edge_ref.first->is_in_domain()) edge_ref = cdt.mirror_edge(edge_ref);
  CGAL_assertion(edge_ref.first->is_in_domain());

  bool flip_ori = ( edge_ref.first->vertex( (edge_ref.second+1)%3 )->info() != source(href, tm) );

  CGAL_assertion(
    ( !flip_ori &&
      edge_ref.first->vertex( (edge_ref.second+1)%3 )->info() == source(href, tm) &&
      edge_ref.first->vertex( (edge_ref.second+2)%3 )->info() == target(href, tm) )
    ||
    ( flip_ori &&
      edge_ref.first->vertex( (edge_ref.second+2)%3 )->info() == source(href, tm) &&
      edge_ref.first->vertex( (edge_ref.second+1)%3 )->info() == target(href, tm) )
  );

  for (auto fh : cdt.finite_face_handles())
  {
    if (!fh->is_in_domain()) continue;

    std::array<vertex_descriptor,3> vrts = make_array(fh->vertex(0)->info(),
                                                      fh->vertex(1)->info(),
                                                      fh->vertex(2)->info());
    if (flip_ori) std::swap(vrts[0], vrts[1]);

    CGAL_assertion(Euler::can_add_face(vrts, tm));
    visitor.before_face_copy(boost::graph_traits<TriangleMesh>::null_face(), tm, tm);
    face_descriptor f = Euler::add_face(vrts, tm);
    visitor.after_face_copy(boost::graph_traits<TriangleMesh>::null_face(), tm, f, tm);
  }
}

template <class Traits,
          class PolygonMesh,
          class VertexPointMap,
          class Visitor>
bool close(PolygonMesh& pm, VertexPointMap vpm, typename Traits::Vector_3 plane_normal, Visitor& visitor)
{
  //using vertex_descriptor = typename boost::graph_traits<PolygonMesh>::vertex_descriptor;
  using halfedge_descriptor = typename boost::graph_traits<PolygonMesh>::halfedge_descriptor;
  // using face_descriptor = typename boost::graph_traits<PolygonMesh>::face_descriptor;
  using Point_3 = typename Traits::Point_3;
  using P_traits = Projection_traits_3<Traits>;

  std::vector< std::vector<Point_3> > polygons;
  std::vector< halfedge_descriptor > border_cycles;
  std::vector< Bbox_3 > bboxes;

  typedef CGAL::dynamic_halfedge_property_t<bool> H_tag;
  typename boost::property_map<PolygonMesh, H_tag>::type
    is_hedge_selected = get(H_tag(), pm, false);

  for (halfedge_descriptor h : halfedges(pm))
  {
    if (is_border(h, pm) && get(is_hedge_selected, h)==false)
    {
      border_cycles.push_back(h);
      polygons.emplace_back();
      bboxes.emplace_back();
      for (halfedge_descriptor he : CGAL::halfedges_around_face(h, pm))
      {
        put(is_hedge_selected, he, true);
        polygons.back().push_back(get(vpm, target(he, pm)));
        bboxes.back()+=polygons.back().back().bbox();
      }
    }
  }

  if (bboxes.empty()) return true;

  Bbox_3 gbox;
  for (const Bbox_3& bb : bboxes)
    gbox+=bb;

  // arrange polygons
  int axis = 0;
  if ((gbox.ymax()-gbox.ymin()) > (gbox.xmax()-gbox.xmin())) axis=1;
  if ((gbox.zmax()-gbox.zmin()) > ((gbox.max)(axis)-(gbox.min)(axis))) axis=2;

  std::vector<std::size_t> poly_ids(polygons.size());
  std::iota(poly_ids.begin(), poly_ids.end(), 0);

  std::sort(poly_ids.begin(), poly_ids.end(),
            [&bboxes, axis](std::size_t i, std::size_t j)
            {
              return (bboxes[i].min)(axis) < (bboxes[j].min)(axis);
            });


  std::vector<std::size_t> simply_connected_faces;
  std::vector<bool> handled(poly_ids.size(), false);

  P_traits ptraits(plane_normal);


  for (std::size_t pid=0; pid<poly_ids.size()-1; ++pid)
  {
    std::size_t curr_poly_id = poly_ids[pid];

    if (handled[curr_poly_id]) continue;

    std::vector<std::size_t> nested;
    for (std::size_t npid=pid+1; npid<poly_ids.size(); ++npid)
    {
      std::size_t next_poly_id = poly_ids[npid];

      if (do_overlap(bboxes[curr_poly_id], bboxes[next_poly_id]))
      {
        //TODO: what about tanjencies?
        //TODO: check orientation for predicate
        if (bounded_side_2(polygons[curr_poly_id].begin(), polygons[curr_poly_id].end(),
                           polygons[next_poly_id].front(), ptraits) == ON_BOUNDED_SIDE)
        {
          nested.push_back(next_poly_id);
          handled[next_poly_id]=true;
        }
      }
      else
        break; // no further overlaps
    }

    handled[curr_poly_id] = true;

    if (nested.empty())
      simply_connected_faces.push_back(curr_poly_id);
  }

  if (!handled[poly_ids.back()])
    simply_connected_faces.push_back(poly_ids.back());

  for (std::size_t id : simply_connected_faces)
  {
    halfedge_descriptor h = border_cycles[id];
    visitor.before_face_copy(boost::graph_traits<PolygonMesh>::null_face(), pm, pm);
    Euler::fill_hole(h, pm);
    visitor.after_face_copy(boost::graph_traits<PolygonMesh>::null_face(), pm, face(h, pm), pm);
  }

  return simply_connected_faces.size()==poly_ids.size();
}

#endif


template <class Plane_3,
          class TriangleMesh,
          class NamedParameters>
Oriented_side
clip_to_bbox(const Plane_3& plane,
             const Bbox_3& bbox,
                   TriangleMesh& tm_out,
             const NamedParameters& np)
{
  typedef typename GetGeomTraits<TriangleMesh, NamedParameters>::type Geom_traits;
  typedef typename Geom_traits::Point_3 Point_3;
  typedef typename GetVertexPointMap<TriangleMesh,
                                     NamedParameters>::type Vpm;

  Vpm vpm_out = parameters::choose_parameter(parameters::get_parameter(np, internal_np::vertex_point),
                                             get_property_map(boost::vertex_point, tm_out));


  std::vector<Point_3> corners(8);
  corners[0] = Point_3(bbox.xmin(),bbox.ymin(),bbox.zmin());
  corners[1] = Point_3(bbox.xmin(),bbox.ymax(),bbox.zmin());
  corners[2] = Point_3(bbox.xmax(),bbox.ymax(),bbox.zmin());
  corners[3] = Point_3(bbox.xmax(),bbox.ymin(),bbox.zmin());
  corners[4] = Point_3(bbox.xmin(),bbox.ymin(),bbox.zmax());
  corners[5] = Point_3(bbox.xmin(),bbox.ymax(),bbox.zmax());
  corners[6] = Point_3(bbox.xmax(),bbox.ymax(),bbox.zmax());
  corners[7] = Point_3(bbox.xmax(),bbox.ymin(),bbox.zmax());

  std::array<CGAL::Oriented_side,8> orientations = {{
    plane.oriented_side(corners[0]),
    plane.oriented_side(corners[1]),
    plane.oriented_side(corners[2]),
    plane.oriented_side(corners[3]),
    plane.oriented_side(corners[4]),
    plane.oriented_side(corners[5]),
    plane.oriented_side(corners[6]),
    plane.oriented_side(corners[7])
  }};

  // description of faces of the bbox
  static constexpr std::array<int, 24> face_indices
    { { 0, 1, 2, 3,
        2, 1, 5, 6,
        3, 2, 6, 7,
        1, 0, 4, 5,
        4, 0, 3, 7,
        6, 5, 4, 7 } };

  static constexpr std::array<int, 24> edge_indices
    { { 0,  1,  2, 3,
        1,  4,  5, 6,
        2,  6,  7, 8,
        0,  9, 10, 4,
        9,  3,  8, 11,
        5, 10, 11, 7 } };

  std::array<int, 12> edge_ipt_id;
  edge_ipt_id.fill(-1);

  auto inter_pt_index =
    [&plane, &corners, &edge_ipt_id](int i, int j, int edge_id)
  {
    if (edge_ipt_id[edge_id]==-1)
    {
      edge_ipt_id[edge_id] = static_cast<int> (corners.size());
      corners.push_back(typename Geom_traits::Construct_plane_line_intersection_point_3()
                      (plane, corners[i], corners[j]));
    }

    return edge_ipt_id[edge_id];
  };

  std::vector< std::vector<int> > output_faces(6);
  bool all_in = true;
  bool all_out = true;
  std::bitset<14> in_point_bits; // to collect the set of points in the clipped bbox

  // for each face of the bbox, we look for intersection of the plane with its edges
  for(int i=0; i<6; ++i)
  {
    for(int k=0; k< 4; ++k)
    {
      int current_id = face_indices[4*i + k];
      int next_id = face_indices[4*i + (k+1)%4];
      int edge_id = edge_indices[4 * i + k];

      switch(orientations[ current_id ])
      {
        case ON_NEGATIVE_SIDE:
        {
          all_out=false;
          // point on or on the negative side
          output_faces[i].push_back(current_id);
          in_point_bits.set(output_faces[i].back());
          // check for intersection of the edge
          if(orientations[ next_id ] == ON_POSITIVE_SIDE)
          {
            output_faces[i].push_back(
              inter_pt_index(current_id, next_id, edge_id));
            in_point_bits.set(output_faces[i].back());
          }
          break;
        }
        case ON_POSITIVE_SIDE:
        {
          all_in = false;
          // check for intersection of the edge
          if(orientations[ next_id ] == ON_NEGATIVE_SIDE)
          {
            output_faces[i].push_back(
              inter_pt_index(current_id, next_id, edge_id));
            in_point_bits.set(output_faces[i].back());
          }
          break;
        }
        case ON_ORIENTED_BOUNDARY:
        {
          output_faces[i].push_back(current_id);
          in_point_bits.set(output_faces[i].back());
        }
      }
    }
    if(output_faces[i].size() < 3){
      CGAL_assertion(output_faces[i].empty() ||
                     (output_faces[i].front()<8 && output_faces[i].back()<8));
      output_faces[i].clear(); // edge of the bbox included in the plane
    }
  }

  // the intersection is the full bbox
  if(all_in) return ON_NEGATIVE_SIDE;
  if(all_out) return ON_POSITIVE_SIDE;

  // build the clipped bbox
  typedef boost::graph_traits<TriangleMesh> graph_traits;
  typedef typename graph_traits::vertex_descriptor vertex_descriptor;
  typedef typename graph_traits::halfedge_descriptor halfedge_descriptor;
  typedef typename graph_traits::face_descriptor face_descriptor;

  std::map<int, vertex_descriptor> out_vertices;
  for(int i=0; i<14;++i)
  {
    if (in_point_bits.test(i))
    {
      vertex_descriptor v = add_vertex(tm_out);
      out_vertices.insert(std::make_pair(i, v));
      put(vpm_out, v, corners[i]);
    }
  }

  std::map< std::pair<int,int>, halfedge_descriptor> hedge_map;
  const halfedge_descriptor null_hedge = graph_traits::null_halfedge();
  const face_descriptor null_fd = graph_traits::null_face();
  for(const std::vector<int>& findices : output_faces)
  {
    if(findices.empty()) continue;
    const face_descriptor fd=add_face(tm_out);
    int prev_id = findices.back();

    // create of recover face boundary halfedges
    std::vector<halfedge_descriptor> hedges;
    hedges.reserve(findices.size());
    for(int current_id : findices)
    {
      vertex_descriptor src = out_vertices[prev_id], tgt = out_vertices[current_id];

      std::pair<typename std::map< std::pair<int,int>,
                halfedge_descriptor>::iterator, bool> res =
        hedge_map.insert(std::make_pair(std::make_pair(prev_id, current_id), null_hedge));
      if(res.second)
      {
        res.first->second = halfedge(add_edge(tm_out), tm_out);
        hedge_map.insert(std::make_pair(std::make_pair(current_id, prev_id),
                            opposite(res.first->second, tm_out)));
        set_face(opposite(res.first->second, tm_out), null_fd, tm_out);

      }
      hedges.push_back(res.first->second);
      // set edge source and target
      set_target(hedges.back(), tgt, tm_out);
      set_target(opposite(hedges.back(), tm_out), src, tm_out);
      // set face pointer of halfedges
      set_face(hedges.back(), fd, tm_out);
      // set vertex halfedge
      set_halfedge(src, opposite(hedges.back(), tm_out), tm_out);
      set_halfedge(tgt, hedges.back(), tm_out);

      if(current_id==findices.front())
        set_halfedge(fd, hedges.back(), tm_out);

      prev_id = current_id;
    }
    CGAL_assertion(hedges.size() == findices.size());

    // set next/prev relationship
    halfedge_descriptor prev_h=hedges.back();
    for(halfedge_descriptor h : hedges)
    {
      set_next(prev_h, h, tm_out);
      prev_h = h;
    }
  }

  // handle the face of the plane:
  // look for a border halfedge and reconstruct the face of the plane
  // by turning around vertices inside the mesh constructed above
  // until we reach another border halfedge
  for(halfedge_descriptor h : halfedges(tm_out))
  {
    if(face(h, tm_out) == null_fd)
    {
      face_descriptor fd = add_face(tm_out);
      set_halfedge(fd, h, tm_out);

      halfedge_descriptor h_prev=h;
      halfedge_descriptor h_curr=h;
      do{
        h_curr=opposite(h_curr, tm_out);
        do{
          h_curr=opposite(prev(h_curr, tm_out), tm_out);
        } while(face(h_curr, tm_out) != null_fd && h_curr!=h);
        set_face(h_prev, fd, tm_out);
        set_next(h_prev, h_curr, tm_out);
        if(h_curr==h)
          break;
        h_prev=h_curr;
      } while(true);
      break;
    }
  }
  CGAL_assertion(is_valid_polygon_mesh(tm_out));

  // triangulate the faces
  CGAL::Polygon_mesh_processing::triangulate_faces(tm_out, np);

  return ON_ORIENTED_BOUNDARY;
}

template <class TriangleMesh, class Ecm, class VPM, class UserVisitor>
void split_along_edges(TriangleMesh& tm,
                       Ecm ecm,
                       VPM vpm,
                       UserVisitor& user_visitor)
{
  typedef boost::graph_traits<TriangleMesh> GT;
  typedef typename GT::face_descriptor face_descriptor;
  typedef typename GT::edge_descriptor edge_descriptor;
  typedef typename GT::vertex_descriptor vertex_descriptor;
  typedef typename GT::halfedge_descriptor halfedge_descriptor;

  std::vector<edge_descriptor> shared_edges;
  for(edge_descriptor e : edges(tm))
    if(get(ecm, e))
      shared_edges.push_back(e);

  std::size_t nb_shared_edges = shared_edges.size();
  std::vector<halfedge_descriptor> hedges_to_update;

  typedef CGAL::dynamic_halfedge_property_t<bool> H_tag;
  typename boost::property_map<TriangleMesh, H_tag>::type
    no_target_update = get(H_tag(), tm);

  std::vector< std::pair<halfedge_descriptor, vertex_descriptor> > vertices_to_duplicate;

  //collect border halfedges having as target one of the edge endpoints
  std::set<halfedge_descriptor> extra_border_hedges;
  for(std::size_t k=0; k<nb_shared_edges; ++k)
  {
    if (is_border(shared_edges[k], tm)) continue;
    for(halfedge_descriptor h : halfedges_around_target(target(shared_edges[k], tm), tm))
      if(is_border(h, tm))
        extra_border_hedges.insert(h);
    for(halfedge_descriptor h : halfedges_around_target(source(shared_edges[k], tm), tm))
      if(is_border(h, tm))
        extra_border_hedges.insert(h);
  }

  for(halfedge_descriptor h : extra_border_hedges)
  {
    put(no_target_update, h, true);
    set_halfedge(target(h, tm), h, tm);
    hedges_to_update.push_back(h);
  }

  // now duplicate the edge and set its pointers
  for(std::size_t k=0; k<nb_shared_edges; ++k)
  {
    if (is_border(shared_edges[k], tm)) continue;
    halfedge_descriptor h    = halfedge(shared_edges[k], tm);
    face_descriptor fh = face(h, tm);
    //add edge
    user_visitor.before_edge_duplicated(h, tm);
    halfedge_descriptor new_hedge = halfedge(add_edge(tm), tm),
                        new_opp   = opposite(new_hedge,tm);
    user_visitor.after_edge_duplicated(h, new_hedge, tm);

    vertex_descriptor vt = target(h, tm);
    vertex_descriptor vs = source(h, tm);

    //replace h with new_hedge
    set_next(new_hedge, next(h, tm), tm);
    set_next(prev(h, tm), new_hedge, tm);
    set_face(new_hedge, fh, tm);
    set_halfedge(fh, new_hedge, tm);

    set_target(new_hedge, vt, tm);
    set_target(new_opp, vs, tm);

    set_face(new_opp, GT::null_face(), tm);
    set_face(h, GT::null_face(), tm);

    // handle vertices to duplicate
    halfedge_descriptor h_vt = halfedge(vt, tm);
    if(get(no_target_update, h_vt))
      vertices_to_duplicate.push_back(std::make_pair(h, vt));
    else
      set_halfedge(vt, h, tm);
    halfedge_descriptor h_vs = halfedge(vs, tm);
    if(get(no_target_update, h_vs))
      vertices_to_duplicate.push_back(std::make_pair(new_opp, vs));
    else
      set_halfedge(vs, new_opp, tm);

    hedges_to_update.push_back(h);
    put(no_target_update, h, true);
    hedges_to_update.push_back(new_opp);
    put(no_target_update, new_opp, true);

    CGAL_assertion(next(prev(new_hedge, tm), tm) == new_hedge);
    CGAL_assertion(prev(next(new_hedge, tm), tm) == new_hedge);
  }

  // update next/prev relationship
  for(halfedge_descriptor h : hedges_to_update)
  {
    CGAL_assertion(is_border(h, tm));
    halfedge_descriptor h_opp = opposite(h, tm);

    // set next pointer of h, visiting faces inside the patch we consider
    halfedge_descriptor candidate = opposite(prev(h_opp, tm), tm);
    while (!is_border(candidate, tm))
      candidate = opposite(prev(candidate, tm), tm);
    set_next(h, candidate, tm);
    CGAL_assertion(prev(next(h_opp, tm), tm)==h_opp);

    CGAL_assertion(prev(next(h, tm), tm) == h);
    CGAL_assertion(is_border(next(h, tm), tm));
  }

  for(const std::pair<halfedge_descriptor, vertex_descriptor>& p : vertices_to_duplicate)
  {
    user_visitor.before_vertex_copy(p.second, tm, tm);
    vertex_descriptor nv = add_vertex(tm);
    put(vpm, nv, get(vpm, p.second));
    user_visitor.after_vertex_copy(p.second, tm, nv, tm);

    for(halfedge_descriptor h : halfedges_around_target(p.first, tm))
      set_target(h, nv, tm);
    set_halfedge(nv, p.first, tm);
   }

  CGAL_assertion_code(for(halfedge_descriptor h : hedges_to_update))
  {
    CGAL_assertion(next(prev(h, tm), tm) == h);
    CGAL_assertion(prev(next(h, tm), tm) == h);
  }

  for(halfedge_descriptor h : hedges_to_update)
  {
    for(halfedge_descriptor hh : halfedges_around_target(h, tm))
        if(h!=hh)
          set_target(hh, target(h, tm), tm);
  }

  CGAL_assertion(is_valid_polygon_mesh(tm));
}

template <class TriangleMesh,
          class NamedParameters1,
          class NamedParameters2>
void
generic_clip_impl(
        TriangleMesh& tm1,
        TriangleMesh& tm2,
  const NamedParameters1& np1,
  const NamedParameters2& np2)
{
  using parameters::choose_parameter;
  using parameters::get_parameter;

// Vertex point maps
  //for input meshes
  typedef typename GetVertexPointMap<TriangleMesh,
                                     NamedParameters1>::type Vpm;
  typedef typename GetVertexPointMap<TriangleMesh,
                                     NamedParameters2>::type Vpm2;

  static_assert(std::is_same<typename boost::property_traits<Vpm>::value_type,
                             typename boost::property_traits<Vpm>::value_type>::value);

  Vpm vpm1 = choose_parameter(get_parameter(np1, internal_np::vertex_point),
                              get_property_map(boost::vertex_point, tm1));

  Vpm2 vpm2 = choose_parameter(get_parameter(np2, internal_np::vertex_point),
                               get_property_map(boost::vertex_point, tm2));

  if (&tm1==&tm2)
  {
    // TODO mark all edges
    return;
  }

  // handle case of empty meshes (isolated vertices are ignored)
  if (faces(tm1).empty())
    return;

// Edge is-constrained maps
  //for input meshes
  typedef typename internal_np::Lookup_named_param_def <
    internal_np::edge_is_constrained_t,
    NamedParameters1,
    Corefinement::No_mark<TriangleMesh>//default
  > ::type User_ecm1;

  // User and internal edge is-constrained map
  typedef typename boost::template property_map<TriangleMesh, CGAL::dynamic_edge_property_t<bool> >::type Algo_ecm1;
  typedef Corefinement::No_mark<TriangleMesh> Ecm2;
  typedef OR_property_map<Algo_ecm1, User_ecm1> Ecm1;
  typedef Corefinement::Ecm_bind<TriangleMesh, Ecm1, Ecm2> Ecm_in;

  Algo_ecm1 algo_ecm1  = get(CGAL::dynamic_edge_property_t<bool>(), tm1);
  Ecm1 ecm1 = Ecm1(algo_ecm1, choose_parameter<User_ecm1>(get_parameter(np1, internal_np::edge_is_constrained)));
  Ecm2 ecm2;

  // Face index point maps
  typedef typename CGAL::GetInitializedFaceIndexMap<TriangleMesh, NamedParameters1>::type FaceIndexMap1;
  FaceIndexMap1 fid_map1 = get_initialized_face_index_map(tm1, np1);

  const bool use_compact_clipper =
    choose_parameter(get_parameter(np1, internal_np::use_compact_clipper), true);


  // User visitor
  typedef typename internal_np::Lookup_named_param_def <
    internal_np::visitor_t,
    NamedParameters1,
    Corefinement::Default_visitor<TriangleMesh>//default
  > ::type User_visitor;
  User_visitor uv(choose_parameter<User_visitor>(get_parameter(np1, internal_np::visitor)));

  // surface intersection algorithm call
  typedef Corefinement::Generic_clip_output_builder<TriangleMesh,
                                                    Vpm, Vpm2,
                                                    Algo_ecm1,
                                                    FaceIndexMap1,
                                                    Default> Ob;

  typedef Corefinement::Surface_intersection_visitor_for_corefinement<
    TriangleMesh, Vpm, Vpm2, Ob, Ecm_in, User_visitor> Algo_visitor;
  Ecm_in ecm_in(tm1,tm2,ecm1,ecm2);
  Ob ob(tm1, tm2, vpm1, vpm2, algo_ecm1, fid_map1, use_compact_clipper);

  Corefinement::Intersection_of_triangle_meshes<TriangleMesh, Vpm, Vpm2, Algo_visitor >
    functor(tm1, tm2, vpm1, vpm2, Algo_visitor(uv,ob,ecm_in,&tm2), &tm2);
  functor(CGAL::Emptyset_iterator(), false, true);
}




#ifndef CGAL_PLANE_CLIP_DO_NOT_USE_TRIANGULATION
template <class PolygonMesh, class Clip_visitor>
struct Visitor_wrapper_for_triangulate_face
  : public ::CGAL::Polygon_mesh_processing::Hole_filling::Default_visitor
//TODO: @afabri --> I shouldn't be the one doing the inheritance...
//TODO: @afabri --> nothing on edge?
{
  using face_descriptor = typename boost::graph_traits<PolygonMesh>::face_descriptor;

  Clip_visitor& clip_visitor;
  const PolygonMesh& pm;
  std::vector<face_descriptor> triangulation_faces;

  Visitor_wrapper_for_triangulate_face(const PolygonMesh& pm, Clip_visitor& clip_visitor)
    : clip_visitor(clip_visitor)
    , pm(pm)
  {}

  void before_subface_creations(face_descriptor f_split)
  {
    CGAL_assertion(triangulation_faces.empty());
    triangulation_faces.push_back(f_split);
  }
  void after_subface_creations()
  {
    for (face_descriptor f : triangulation_faces)
    {
      clip_visitor.before_face_copy(boost::graph_traits<PolygonMesh>::null_face(), pm, pm);
      clip_visitor.after_face_copy(boost::graph_traits<PolygonMesh>::null_face(), pm, f, pm);
    }
    triangulation_faces.clear();
  }

  void after_subface_created(face_descriptor f_new)
  {
    triangulation_faces.push_back(f_new);
  }
};
#endif

} // end of internal namespace

/**
  * \ingroup PMP_corefinement_grp
  *
  * \brief clips `tm` by keeping the part that is inside the volume \link coref_def_subsec bounded \endlink by `clipper`.
  *
  * If `tm` is closed, the clipped part can be kept closed by setting the named parameter `clip_volume` to `true`.
  * See Subsection \ref coref_clip for more details.
  *
  * \attention With the current implementation, `clipper` will be modified (refined with the intersection with `tm`).
  *
  * \pre \link CGAL::Polygon_mesh_processing::does_self_intersect() `!CGAL::Polygon_mesh_processing::does_self_intersect(tm)` \endlink
  * \pre \link CGAL::Polygon_mesh_processing::does_self_intersect() `!CGAL::Polygon_mesh_processing::does_self_intersect(clipper)` \endlink
  * \pre \link CGAL::Polygon_mesh_processing::does_bound_a_volume() `CGAL::Polygon_mesh_processing::does_bound_a_volume(clipper)` \endlink
  *
  * @tparam TriangleMesh a model of `MutableFaceGraph`, `HalfedgeListGraph` and `FaceListGraph`.
  * @tparam NamedParameters1 a sequence of \ref bgl_namedparameters "Named Parameters"
  * @tparam NamedParameters2 a sequence of \ref bgl_namedparameters "Named Parameters"
  *
  * @param tm input triangulated surface mesh
  * @param clipper triangulated surface mesh used to clip `tm`
  * @param np_tm an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
  * @param np_c an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
  *
  *
  * \cgalNamedParamsBegin
  *   \cgalParamNBegin{vertex_point_map}
  *     \cgalParamDescription{a property map associating points to the vertices of `tm` (resp. `clipper`)}
  *     \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
  *                    as key type and `%Point_3` as value type}
  *     \cgalParamDefault{`boost::get(CGAL::vertex_point, tm (resp. clipper))`}
  *     \cgalParamExtra{If this parameter is omitted, an internal property map for `CGAL::vertex_point_t`
  *                     must be available in `TriangleMesh`.}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{face_index_map}
  *     \cgalParamDescription{a property map associating to each face of `tm` (`clipper`) a unique index between `0` and `num_faces(tm (resp. clipper)) - 1`}
  *     \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<TriangleMesh>::%face_descriptor`
  *                    as key type and `std::size_t` as value type}
  *     \cgalParamDefault{an automatically indexed internal map}
  *     \cgalParamExtra{if the property map is writable, the indices of the faces of `tm` and `clipper`
  *                     will be set after refining `tm` with the intersection with `clipper`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{visitor}
  *     \cgalParamDescription{a visitor used to track the creation of new faces}
  *     \cgalParamType{a class model of `PMPCorefinementVisitor`}
  *     \cgalParamDefault{`Corefinement::Default_visitor<TriangleMesh>`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{throw_on_self_intersection}
  *     \cgalParamDescription{If `true`, the set of triangles close to the intersection of `tm` and `clipper` will be
  *                           checked for self-intersections and `Corefinement::Self_intersection_exception`
  *                           will be thrown if at least one self-intersection is found.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{clip_volume}
  *     \cgalParamDescription{If `true`, and `tm` is closed, the clipping will be done on the volume
  *                           \link coref_def_subsec bounded \endlink by `tm` rather than on its surface
  *                           (i.e., `tm` will be kept closed).}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{use_compact_clipper}
  *     \cgalParamDescription{if `false`, the parts of `tm` coplanar with `clipper` will not be part of the output.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`true`}
  *     \cgalParamExtra{This option has an effect only if a surface and not a volume is clipped,
  *                     (i.e., if `clip_volume` is `false` or if `tm` is open).}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{do_not_modify}
  *     \cgalParamDescription{(`np_c` only) if `true`, `clipper` will not be modified.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *     \cgalParamExtra{If this option is set to `true`, `tm` is no longer required to be without self-intersection.
  *                     Setting this option to `true` will automatically set `throw_on_self_intersection` to `false`
  *                     and `clip_volume` to `false`.}
  *   \cgalParamNEnd
  * \cgalNamedParamsEnd
  *
  * @return `true` if the output surface mesh is manifold.
  *         If `false` is returned `tm` and `clipper` are only corefined.
  *
  * @see `split()`
  */
template <class TriangleMesh,
          class NamedParameters1 = parameters::Default_named_parameters,
          class NamedParameters2 = parameters::Default_named_parameters>
bool
clip(TriangleMesh& tm,
     TriangleMesh& clipper,
     const NamedParameters1& np_tm = parameters::default_values(),
     const NamedParameters2& np_c = parameters::default_values())
{
  if (parameters::choose_parameter(parameters::get_parameter(np_c, internal_np::do_not_modify), false))
  {
    CGAL_assertion(is_closed(clipper));

    internal::generic_clip_impl(tm, clipper, np_tm, np_c);
    return true;
  }

  const bool clip_volume =
    parameters::choose_parameter(parameters::get_parameter(np_tm, internal_np::clip_volume), false);

  if(clip_volume && is_closed(tm))
    return corefine_and_compute_intersection(tm, clipper, tm, np_tm, np_c);
  return corefine_and_compute_intersection(tm, clipper, tm,
                                           np_tm.use_bool_op_to_clip_surface(true),
                                           np_c);
}

/**
  * \ingroup PMP_corefinement_grp
  *
  * \brief clips `pm` by keeping the part that is on the negative side of `plane` (the side opposite to its normal vector).
  *
  * If `pm` is closed, the clipped part can be kept closed by setting the named parameter `clip_volume`to `true`.
  * See Subsection \ref coref_clip for more details.
  *
  * @tparam PolygonMesh a model of `MutableFaceGraph`, `HalfedgeListGraph` and `FaceListGraph`.
  *                      An internal property map for `CGAL::vertex_point_t` must be available.
  *
  * @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
  *
  * @param pm input surface mesh
  * @param plane plane whose negative side defines the halfspace to intersect `pm` with.
  *              `Plane_3` is the plane type for the same CGAL kernel as the point of the vertex point map of `pm`.
  * @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
  *
  * \cgalNamedParamsBegin
  *
  *   \cgalParamNBegin{concurrency_tag}
  *     \cgalParamDescription{a tag indicating if the task should be performed using one or several threads.}
  *     \cgalParamType{Either `CGAL::Sequential_tag`, or `CGAL::Parallel_tag`, or `CGAL::Parallel_if_available_tag`}
  *     \cgalParamDefault{`CGAL::Sequential_tag`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{vertex_point_map}
  *     \cgalParamDescription{a property map associating points to the vertices of `pm`}
  *     \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<PolygonMesh>::%vertex_descriptor`
  *                    as key type and `%Point_3` as value type}
  *     \cgalParamDefault{`boost::get(CGAL::vertex_point, pm)`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{visitor}
  *     \cgalParamDescription{a visitor used to track the creation of new faces, edges, and faces.
  *                           Note that as there are no mesh associated with `plane`,
  *                           `boost::graph_traits<PolygonMesh>::null_halfedge()` and `boost::graph_traits<PolygonMesh>::null_face()` will be used when calling
  *                           functions of the visitor expecting a halfedge or a face from `plane`. Similarly, `pm` will be used as the mesh of `plane`.}
  *     \cgalParamType{a class model of `PMPCorefinementVisitor`}
  *     \cgalParamDefault{`Corefinement::Default_visitor<PolygonMesh>`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{throw_on_self_intersection}
  *     \cgalParamDescription{If `true`, the set of triangles close to the intersection of `pm`
  *                           and `plane` will be checked for self-intersections
  *                           and `CGAL::Polygon_mesh_processing::Corefinement::Self_intersection_exception`
  *                           will be thrown if at least one self-intersection is found.
  *                           This option is only taken into account if `pm` is a triangle mesh.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{allow_self_intersections}
  *     \cgalParamDescription{If `true`, self-intersections in `pm` are accepted.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *     \cgalParamExtra{If this option is set to `true`, `pm` is no longer required to be without self-intersection.
  *                     Setting this option to `true` will automatically set `throw_on_self_intersection` to `false`
  *                     and `clip_volume` to `false` (overwriting any value provided)}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{clip_volume}
  *     \cgalParamDescription{If `true`, and if `pm` is closed, the clipping will be done on
  *                           the volume \link coref_def_subsec bounded \endlink by `pm`
  *                           rather than on its surface (i.e., `pm` will remain closed).}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{use_compact_clipper}
  *     \cgalParamDescription{If `false`, the parts of `pm` coplanar with `plane` will not be part of the output.
  *                           Always `true` if `clip_volume` is `true`.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`true`}
  *   \cgalParamNEnd
  *
  *    \cgalParamNBegin{do_not_triangulate_faces}
  *      \cgalParamDescription{If the input mesh is triangulated and this parameter is set to `false`, the mesh will be kept triangulated.
  *                            Always `true` if `pm` is not a triangle mesh.}
  *      \cgalParamType{Boolean}
  *      \cgalParamDefault{`false`}
  *    \cgalParamNEnd
  *
  * \cgalNamedParamsEnd
  *
  * @return `true`
  *
  * @see `split()`
  */
template <class PolygonMesh,
          class NamedParameters = parameters::Default_named_parameters>
bool clip(PolygonMesh& pm,
#ifdef DOXYGEN_RUNNING
          const Plane_3& plane,
#else
          const typename GetGeomTraits<PolygonMesh, NamedParameters>::type::Plane_3& plane,
#endif
          const NamedParameters& np = parameters::default_values())
{
  using parameters::choose_parameter;
  using parameters::get_parameter;
  using parameters::get_parameter_reference     ;

  using halfedge_descriptor = typename boost::graph_traits<PolygonMesh>::halfedge_descriptor;

  using GT = typename GetGeomTraits<PolygonMesh, NamedParameters>::type;
  GT traits = choose_parameter<GT>(get_parameter(np, internal_np::geom_traits));
  auto vpm = choose_parameter(get_parameter(np, internal_np::vertex_point),
                              get_property_map(vertex_point, pm));


  using Default_visitor = Corefinement::Default_visitor<PolygonMesh>;
  Default_visitor default_visitor;
  using Visitor_ref = typename internal_np::Lookup_named_param_def<internal_np::visitor_t, NamedParameters, Default_visitor>::reference;
  Visitor_ref visitor = choose_parameter(get_parameter_reference(np, internal_np::visitor), default_visitor);
  // constexpr bool has_visitor = !std::is_same_v<Default_visitor, std::remove_cv_t<std::remove_reference_t<Visitor_ref>>>;

  typedef typename internal_np::Lookup_named_param_def <
    internal_np::concurrency_tag_t,
    NamedParameters,
    Sequential_tag
  > ::type Concurrency_tag;

  // config flags
  bool clip_volume =
    parameters::choose_parameter(parameters::get_parameter(np, internal_np::clip_volume), false);
  bool use_compact_clipper =
    choose_parameter(get_parameter(np, internal_np::use_compact_clipper), true);
  const bool throw_on_self_intersection =
    choose_parameter(get_parameter(np, internal_np::throw_on_self_intersection), false);
  const bool allow_self_intersections =
    choose_parameter(get_parameter(np, internal_np::allow_self_intersections), false);
  bool triangulate = !choose_parameter(get_parameter(np, internal_np::do_not_triangulate_faces), false);

  auto vos = get(dynamic_vertex_property_t<Oriented_side>(), pm);
  auto ecm = get(dynamic_edge_property_t<bool>(), pm, false);

  if (triangulate && !is_triangle_mesh(pm))
    triangulate = false;

  refine_with_plane(pm, plane, parameters::vertex_oriented_side_map(vos)
                                          .edge_is_marked_map(ecm)
                                          .vertex_point_map(vpm)
                                          .geom_traits(traits)
                                          .do_not_triangulate_faces(!triangulate)
                                          .throw_on_self_intersection(!allow_self_intersections &&
                                                                      throw_on_self_intersection)
                                          .visitor(std::ref(visitor))
                                          .concurrency_tag(Concurrency_tag()));

  CGAL_assertion(is_valid_polygon_mesh(pm));

  if (allow_self_intersections)
    clip_volume=false;

  if (clip_volume && !is_closed(pm)) clip_volume=false;
  if (clip_volume && use_compact_clipper) use_compact_clipper=false;

  CGAL_precondition(!clip_volume || !triangulate || does_bound_a_volume(pm));

  auto fcc = get(dynamic_face_property_t<std::size_t>(), pm);

  std::size_t nbcc = connected_components(pm, fcc, CGAL::parameters::edge_is_constrained_map(ecm));

  std::vector<bool> classified(nbcc, false);
  std::vector<std::size_t> ccs_to_remove;

  for (auto f : faces(pm))
  {
    std::size_t ccid = get(fcc, f);
    if (classified[ccid]) continue;
    halfedge_descriptor hf = halfedge(f, pm);
    for(halfedge_descriptor h : CGAL::halfedges_around_face(hf, pm))
    {
      CGAL::Oriented_side os = get(vos, target(h,pm));
      if (os==CGAL::ON_ORIENTED_BOUNDARY) continue;
      classified[ccid]=true;
      if (os==CGAL::ON_POSITIVE_SIDE) ccs_to_remove.push_back(ccid);
    }

    if (!classified[ccid])
    {
      if (!use_compact_clipper) ccs_to_remove.push_back(ccid);
      classified[ccid]=true;
    }
  }

  remove_connected_components(pm, ccs_to_remove, fcc);

  if (clip_volume)
  {
    //TODO: add in the traits construct_orthogonal_vector
    if (triangulate)
      internal::close_and_triangulate<GT>(pm, vpm, plane.orthogonal_vector(), visitor);
    else
      if (!internal::close<GT>(pm, vpm, plane.orthogonal_vector(), visitor))
        internal::close_and_triangulate<GT>(pm, vpm, plane.orthogonal_vector(), visitor);
  }

  return true;
}

/**
  * \ingroup PMP_corefinement_grp
  *
  * \brief clips `tm` by keeping the part that is inside `iso_cuboid`.
  *
  * If `tm` is closed, the clipped part can be kept closed by setting the named parameter `clip_volume` to `true`.
  * See Subsection \ref coref_clip for more details.
  *
  * \note `Iso_cuboid_3` must be from the same kernel as the point of the internal vertex point map of `TriangleMesh`.
  * \note `Iso_cuboid_3` must be from the same kernel as the point of the vertex point map of `tm`.
  *
  * \pre \link CGAL::Polygon_mesh_processing::does_self_intersect() `!CGAL::Polygon_mesh_processing::does_self_intersect(tm)` \endlink
  *
  * @tparam TriangleMesh a model of `MutableFaceGraph`, `HalfedgeListGraph` and `FaceListGraph`.
  *                      An internal property map for `CGAL::vertex_point_t` must be available.
  *
  * @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
  *
  * @param tm input triangulated surface mesh
  * @param iso_cuboid iso-cuboid used to clip `tm`.
  * @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
  *
  * \cgalNamedParamsBegin
  *   \cgalParamNBegin{visitor}
  *     \cgalParamDescription{a visitor used to track the creation of new faces}
  *     \cgalParamType{a class model of `PMPCorefinementVisitor`}
  *     \cgalParamDefault{`Corefinement::Default_visitor<TriangleMesh>`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{clip_volume}
  *     \cgalParamDescription{If `true`, and `tm` is closed, the clipping will be done on
  *                           the volume \link coref_def_subsec bounded \endlink by `tm`
  *                           rather than on its surface (i.e., `tm` will be kept closed).}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{use_compact_clipper}
  *     \cgalParamDescription{if `false` the parts of `tm` coplanar with `iso_cuboid` will not be part of the output}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`true`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{throw_on_self_intersection}
  *     \cgalParamDescription{If `true`, the set of triangles close to the intersection of `tm`
  *                           and `iso_cuboid` will be checked for self-intersections
  *                           and `CGAL::Polygon_mesh_processing::Corefinement::Self_intersection_exception`
  *                           will be thrown if at least one self-intersection is found.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{allow_self_intersections}
  *     \cgalParamDescription{If `true`, self-intersections in `tm` are accepted.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *     \cgalParamExtra{If this option is set to `true`, `tm` is no longer required to be without self-intersection.
  *                     Setting this option to `true` will automatically set `throw_on_self_intersection` to `false`
  *                     and `clip_volume` to `false`.}
  *   \cgalParamNEnd
  *
  * \cgalNamedParamsEnd
  *
  * @return `true` if the output surface mesh is manifold.
  *         If `false` is returned `tm` is only refined by the intersection with `iso_cuboid`.
  *
  * @see `split()`
  */
template <class TriangleMesh,
          class NamedParameters = parameters::Default_named_parameters>
bool clip(TriangleMesh& tm,
#ifdef DOXYGEN_RUNNING
          const Iso_cuboid_3& iso_cuboid,
#else
          const typename GetGeomTraits<TriangleMesh, NamedParameters>::type::Iso_cuboid_3& iso_cuboid,
#endif
          const NamedParameters& np = parameters::default_values())
{
  namespace PMP = CGAL::Polygon_mesh_processing;
  namespace params = CGAL::parameters;

  using params::get_parameter;
  using params::choose_parameter;

  if(std::begin(faces(tm))==std::end(faces(tm))) return true;
  TriangleMesh clipper;

  make_hexahedron(iso_cuboid[0], iso_cuboid[1], iso_cuboid[2], iso_cuboid[3],
                  iso_cuboid[4], iso_cuboid[5], iso_cuboid[6], iso_cuboid[7],
                  clipper);
  triangulate_faces(clipper);

  const bool do_not_modify = choose_parameter(get_parameter(np, internal_np::allow_self_intersections), false);
  return clip(tm, clipper, np, params::do_not_modify(do_not_modify));
}

/*!
  * \ingroup PMP_corefinement_grp
  *
  * corefines `tm` and `splitter` and duplicates edges in `tm` that are on the intersection with `splitter`.
  *
  * \pre \link CGAL::Polygon_mesh_processing::does_self_intersect() `!CGAL::Polygon_mesh_processing::does_self_intersect(tm)` \endlink
  * \pre \link CGAL::Polygon_mesh_processing::does_self_intersect() `!CGAL::Polygon_mesh_processing::does_self_intersect(splitter)` \endlink
  *
  * @tparam TriangleMesh a model of `MutableFaceGraph`, `HalfedgeListGraph`, and `FaceListGraph`.
  *
  * @tparam NamedParameters1 a sequence of \ref bgl_namedparameters "Named Parameters"
  * @tparam NamedParameters2 a sequence of \ref bgl_namedparameters "Named Parameters"
  *
  * @param tm input triangulated surface mesh
  * @param splitter triangulated surface mesh used to split `tm`
  * @param np_tm an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
  * @param np_s an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
  *
  * \cgalNamedParamsBegin
  *   \cgalParamNBegin{vertex_point_map}
  *     \cgalParamDescription{a property map associating points to the vertices of `tm` (`splitter`)}
  *     \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
  *                    as key type and `%Point_3` as value type}
  *     \cgalParamDefault{`boost::get(CGAL::vertex_point, tm)`}
  *     \cgalParamDefault{If this parameter is omitted, an internal property map for
  *                       `CGAL::vertex_point_t` must be available in `TriangleMesh`.}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{visitor}
  *     \cgalParamDescription{(`np_tm` only) a visitor used to track a series of events such as the creation
  *                           of new faces, edges... in both `tm` and `splitter`.}
  *     \cgalParamType{a class model of `PMPCorefinementVisitor`}
  *     \cgalParamDefault{`Corefinement::Default_visitor<TriangleMesh>`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{throw_on_self_intersection}
  *     \cgalParamDescription{If `true`, the set of triangles closed to the intersection of `tm`
  *                           and `splitter` will be checked for self-intersections
  *                           and `CGAL::Polygon_mesh_processing::Corefinement::Self_intersection_exception`
  *                           will be thrown if at least one self-intersection is found.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{do_not_modify}
  *     \cgalParamDescription{(`np_s` only) if `true`, `splitter` will not be modified.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *     \cgalParamExtra{If this option is set to `true`, `tm` is no longer required to be without self-intersection.
  *                     Setting this option to `true` will automatically set `throw_on_self_intersection` to `false`
  *                     and `clip_volume` to `false`.}
  *   \cgalParamNEnd
  *
  * \cgalNamedParamsEnd
  *
  * @see `clip()`
*/
template <class TriangleMesh,
          class NamedParameters1 = parameters::Default_named_parameters,
          class NamedParameters2 = parameters::Default_named_parameters>
void split(TriangleMesh& tm,
           TriangleMesh& splitter,
           const NamedParameters1& np_tm = parameters::default_values(),
           const NamedParameters2& np_s = parameters::default_values())
{
  namespace PMP = CGAL::Polygon_mesh_processing;

  using parameters::get_parameter;
  using parameters::choose_parameter;

  typedef typename GetVertexPointMap<TriangleMesh, NamedParameters1>::type VPM1;
  typedef typename GetVertexPointMap<TriangleMesh, NamedParameters2>::type VPM2;

  typedef typename boost::template property_map<TriangleMesh, CGAL::dynamic_edge_property_t<bool> >::type Ecm;

  VPM1 vpm_tm = choose_parameter(get_parameter(np_tm, internal_np::vertex_point),
                                 get_property_map(vertex_point, tm));
  VPM2 vpm_s = choose_parameter(get_parameter(np_s, internal_np::vertex_point),
                                get_property_map(vertex_point, splitter));

  Ecm ecm  = get(CGAL::dynamic_edge_property_t<bool>(), tm);

  typedef typename internal_np::Lookup_named_param_def <
    internal_np::visitor_t,
    NamedParameters1,
    Corefinement::Default_visitor<TriangleMesh>//default
  > ::type User_visitor;
  User_visitor uv(choose_parameter<User_visitor>(get_parameter(np_tm, internal_np::visitor)));


  // create a constrained edge map and corefine input mesh with the splitter,
  // and mark edges

  const bool do_not_modify_splitter = choose_parameter(get_parameter(np_s, internal_np::do_not_modify), false);

  PMP::corefine(tm, splitter,
                CGAL::parameters::vertex_point_map(vpm_tm).edge_is_constrained_map(ecm).visitor(uv),
                CGAL::parameters::vertex_point_map(vpm_s).do_not_modify(do_not_modify_splitter));

  //split mesh along marked edges
  internal::split_along_edges(tm, ecm, vpm_tm, uv);
}

/**
  * \ingroup PMP_corefinement_grp
  *
  * splits a polygon mesh with a plane.
  *
  * The polygon mesh is refined with the intersection edges, and those edges are duplicated as to create a boundary,
  * and thus separate connected components on either side of the plane.
  *
  * @tparam PolygonMesh a model of `MutableFaceGraph`, `HalfedgeListGraph`, and `FaceListGraph`.
  *                      An internal property map for `CGAL::vertex_point_t` must be available.
  * @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
  *
  * @param pm input surface mesh
  * @param plane the plane that will be used to split `pm`.
  *              `Plane_3` is the plane type for the same CGAL kernel as the point of the vertex point map of `pm`.
  * @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
  *
  * \cgalNamedParamsBegin
  *   \cgalParamNBegin{vertex_point_map}
  *     \cgalParamDescription{a property map associating points to the vertices of `pm`}
  *     \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
  *                    as key type and `%Point_3` as value type}
  *     \cgalParamDefault{`boost::get(CGAL::vertex_point, pm)`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{visitor}
  *     \cgalParamDescription{a visitor used to track the creation of new faces, edges, and vertices.
  *                           Note that as there are no mesh associated with `plane`,
  *                           `boost::graph_traits<PolygonMesh>::null_halfedge()` and `boost::graph_traits<PolygonMesh>::null_face()` will be used when calling
  *                           functions of the visitor expecting a halfedge or a face from `plane`. Similarly, `pm` will be used as the mesh of `plane`.}}
  *     \cgalParamType{a class model of `PMPCorefinementVisitor`}
  *     \cgalParamDefault{`Corefinement::Default_visitor<TriangleMesh>`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{throw_on_self_intersection}
  *     \cgalParamDescription{If `true`, the set of triangles close to the intersection of `pm`
  *                           and `plane` will be checked for self-intersections
  *                           and `CGAL::Polygon_mesh_processing::Corefinement::Self_intersection_exception`
  *                           will be thrown if at least one self-intersection is found.
  *                           This option is only taken into account if `pm` is a triangle mesh.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *    \cgalParamNBegin{do_not_triangulate_faces}
  *      \cgalParamDescription{If the input mesh is triangulated and this parameter is set to `false`, the mesh will be kept triangulated.
  *                            Always `true` if `pm` is not a triangle mesh.}
  *      \cgalParamType{Boolean}
  *      \cgalParamDefault{`false`}
  *    \cgalParamNEnd
  * \cgalNamedParamsEnd
  *
  * @see `clip()`
  */
template <class PolygonMesh,
          class NamedParameters = parameters::Default_named_parameters>
void split(PolygonMesh& pm,
#ifdef DOXYGEN_RUNNING
          const Plane_3& plane,
#else
          const typename GetGeomTraits<PolygonMesh, NamedParameters>::type::Plane_3& plane,
#endif
          const NamedParameters& np = parameters::default_values())
{
  using parameters::choose_parameter;
  using parameters::get_parameter;
  using parameters::get_parameter_reference;

  using GT = typename GetGeomTraits<PolygonMesh, NamedParameters>::type;
  GT traits = choose_parameter<GT>(get_parameter(np, internal_np::geom_traits));
  auto vpm = choose_parameter(get_parameter(np, internal_np::vertex_point),
                              get_property_map(vertex_point, pm));

  typedef typename internal_np::Lookup_named_param_def <
    internal_np::concurrency_tag_t,
    NamedParameters,
    Sequential_tag
  > ::type Concurrency_tag;

  // config flags
  const bool throw_on_self_intersection =
    choose_parameter(get_parameter(np, internal_np::throw_on_self_intersection), false);
  bool triangulate = !choose_parameter(get_parameter(np, internal_np::do_not_triangulate_faces), false);

  auto vos = get(dynamic_vertex_property_t<Oriented_side>(), pm);
  auto ecm = get(dynamic_edge_property_t<bool>(), pm, false);

  if (triangulate && !is_triangle_mesh(pm))
    triangulate = false;

  using Default_visitor = Corefinement::Default_visitor<PolygonMesh>;
  Default_visitor default_visitor;
  using Visitor_ref = typename internal_np::Lookup_named_param_def<internal_np::visitor_t, NamedParameters, Default_visitor>::reference;
  Visitor_ref visitor = choose_parameter(get_parameter_reference(np, internal_np::visitor), default_visitor);

  refine_with_plane(pm, plane, parameters::vertex_oriented_side_map(vos)
                                          .edge_is_marked_map(ecm)
                                          .vertex_point_map(vpm)
                                          .geom_traits(traits)
                                          .do_not_triangulate_faces(!triangulate)
                                          .throw_on_self_intersection(throw_on_self_intersection)
                                          .concurrency_tag(Concurrency_tag())
                                          .visitor(std::ref(visitor)));

  //split mesh along marked edges
  internal::split_along_edges(pm, ecm, vpm, visitor);
}


/**
  * \ingroup PMP_corefinement_grp
  *
  * adds intersection edges of `iso_cuboid` and `tm` in `tm` and duplicates those edges.
  *
  * \note `Iso_cuboid_3` must be from the same kernel as the point of the internal vertex point map of `TriangleMesh`.
  * \note `Iso_cuboid_3` must be from the same kernel as the point of the vertex point map of `tm`.
  *
  * \pre \link CGAL::Polygon_mesh_processing::does_self_intersect() `!CGAL::Polygon_mesh_processing::does_self_intersect(tm)` \endlink
  *
  * @tparam TriangleMesh a model of `MutableFaceGraph`, `HalfedgeListGraph`, and `FaceListGraph`.
  *                      An internal property map for `CGAL::vertex_point_t` must be available.
  * @tparam NamedParameters a sequence of \ref bgl_namedparameters "Named Parameters"
  *
  * @param tm input triangulated surface mesh
  * @param iso_cuboid iso-cuboid used to split `tm`.
  * @param np an optional sequence of \ref bgl_namedparameters "Named Parameters" among the ones listed below
  *
  * \cgalNamedParamsBegin
  *   \cgalParamNBegin{vertex_point_map}
  *     \cgalParamDescription{a property map associating points to the vertices of `tm`}
  *     \cgalParamType{a class model of `ReadWritePropertyMap` with `boost::graph_traits<TriangleMesh>::%vertex_descriptor`
  *                    as key type and `%Point_3` as value type}
  *     \cgalParamDefault{`boost::get(CGAL::vertex_point, tm)`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{visitor}
  *     \cgalParamDescription{a visitor used to track the creation of new faces}
  *     \cgalParamType{a class model of `PMPCorefinementVisitor`}
  *     \cgalParamDefault{`Corefinement::Default_visitor<TriangleMesh>`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{clip_volume}
  *     \cgalParamDescription{If `true`, and `tm` is closed, the clipping will be done on
  *                           the volume \link coref_def_subsec bounded \endlink by `tm`
  *                           rather than on its surface (i.e., `tm` will be kept closed).}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{use_compact_clipper}
  *     \cgalParamDescription{if `false`, the parts of `tm` coplanar with `iso_cuboid` will not be part of the output.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`true`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{throw_on_self_intersection}
  *     \cgalParamDescription{If `true`, the set of triangles close to the intersection of `tm`
  *                           and `iso_cuboid` will be checked for self-intersections
  *                           and `CGAL::Polygon_mesh_processing::Corefinement::Self_intersection_exception`
  *                           will be thrown if at least one self-intersection is found.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *   \cgalParamNEnd
  *
  *   \cgalParamNBegin{allow_self_intersections}
  *     \cgalParamDescription{If `true`, self-intersections in `tm` are accepted.}
  *     \cgalParamType{Boolean}
  *     \cgalParamDefault{`false`}
  *     \cgalParamExtra{If this option is set to `true`, `tm` is no longer required to be without self-intersection.
  *                     Setting this option to `true` will automatically set `throw_on_self_intersection` to `false`
  *                     and `clip_volume` to `false`.}
  *   \cgalParamNEnd
  *
  * \cgalNamedParamsEnd
  *
  * @see `clip()`
  */
template <class TriangleMesh,
          class NamedParameters = parameters::Default_named_parameters>
void split(TriangleMesh& tm,
#ifdef DOXYGEN_RUNNING
           const Iso_cuboid_3& iso_cuboid,
#else
           const typename GetGeomTraits<TriangleMesh, NamedParameters>::type::Iso_cuboid_3& iso_cuboid,
#endif
           const NamedParameters& np = parameters::default_values())
{
  namespace PMP = CGAL::Polygon_mesh_processing;
  namespace params = CGAL::parameters;

  using params::get_parameter;
  using params::choose_parameter;

  TriangleMesh splitter;

  make_hexahedron(iso_cuboid[0], iso_cuboid[1], iso_cuboid[2], iso_cuboid[3],
                  iso_cuboid[4], iso_cuboid[5], iso_cuboid[6], iso_cuboid[7],
                  splitter);
  triangulate_faces(splitter);
  const bool do_not_modify = choose_parameter(get_parameter(np, internal_np::allow_self_intersections), false);
  return split(tm, splitter, np, params::do_not_modify(do_not_modify));
}

} } //end of namespace CGAL::Polygon_mesh_processing

#endif // CGAL_POLYGON_MESH_PROCESSING_CLIP_H