File: createCoupleMatches.C

package info (click to toggle)
openfoam 4.1%2Bdfsg1-1
  • links: PTS, VCS
  • area: main
  • in suites: stretch
  • size: 163,028 kB
  • ctags: 58,990
  • sloc: cpp: 830,760; sh: 10,227; ansic: 8,215; xml: 745; lex: 437; awk: 194; sed: 91; makefile: 77; python: 18
file content (1508 lines) | stat: -rw-r--r-- 58,369 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
1503
1504
1505
1506
1507
1508
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
     \\/     M anipulation  |
-------------------------------------------------------------------------------
License
    This file is part of OpenFOAM.

    OpenFOAM is free software: you can redistribute it and/or modify it
    under the terms of the GNU General Public License as published by
    the Free Software Foundation, either version 3 of the License, or
    (at your option) any later version.

    OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
    ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
    FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
    for more details.

    You should have received a copy of the GNU General Public License
    along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.

Description
    Create coupled match faces and add them to the cells

\*---------------------------------------------------------------------------*/

#include "starMesh.H"
#include "boolList.H"
#include "pointHit.H"
#include "IOmanip.H"
#include "boundBox.H"
#include "Map.H"
#include "unitConversion.H"

// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //

void Foam::starMesh::createCoupleMatches()
{
    // Loop through all couples and create intersection faces. Add all points
    // of intersection faces to the couple points lists. The numbering of
    // the list is set such that the list can be appended to the
    // existing points list

    // Estimate the number of cells affected by couple matches
    const label cellMapSize = min
    (
        cellShapes_.size()/10,
        couples_.size()*2
    );

    // Store newly created faces for each cell
    Map<SLList<face>> cellAddedFaces(cellMapSize);

    Map<SLList<label>> cellRemovedFaces(cellMapSize);

    // In order to remove often allocation, remember the number of live points.
    // If you run out of space in point creation, increase it by the number of
    // couples (good scale) and resize at the end;
    label nLivePoints = points_.size();

    const label infoJump = max(1000, couples_.size()/20);

    forAll(couples_, coupleI)
    {
        if (coupleI % infoJump == 0)
        {
            Info<< "Doing couple " << coupleI << ". STAR couple ID: "
                << couples_[coupleI].coupleID() << endl;
        }

        // Initialise cell edges for master and slave cells
        const coupledFacePair& fp = couples_[coupleI];
        const face& masterFace = cellFaces_[fp.masterCell()][fp.masterFace()];
        const face& slaveFace = cellFaces_[fp.slaveCell()][fp.slaveFace()];

        #ifdef DEBUG_COUPLE
        Info<< "coupleI: " << coupleI << endl
            << "masterFace: " << masterFace << endl
            << "master points: " << masterFace.points(points_) << endl
            << "slaveFace: " << slaveFace << endl
            << "slave points: " << slaveFace.points(points_)
            << endl << endl;
        #endif

        // check the angle of face area vectors
         scalar faceAreaAngle =
             mag
             (
                 -(masterFace.normal(points_) & slaveFace.normal(points_))/
                 (masterFace.mag(points_)*slaveFace.mag(points_) + VSMALL)
             );

         if (faceAreaAngle < 0.94)
         {
             Info<< "Couple direction mismatch in the couple match "
                 << coupleI << ". STAR couple ID: "
                 << couples_[coupleI].coupleID() << endl
                 << "The angle between face normals is "
                 << radToDeg(Foam::acos(faceAreaAngle))
                 << " deg." << endl
                 << "master cell: " << fp.masterCell()
                 << " STAR number: " << starCellID_[fp.masterCell()]
                 << " type: " << cellShapes_[fp.masterCell()].model().name()
                 << " face: " << fp.masterFace() << endl
                 << "slave cell : " << fp.slaveCell()
                 << " STAR number: " << starCellID_[fp.slaveCell()]
                 << " type: " << cellShapes_[fp.slaveCell()].model().name()
                 << " face: " << fp.slaveFace() << endl;
         }

        // Deal with integral patches
        if (fp.integralMatch())
        {
            // Master face is replaced by a set of slave faces

            Map<SLList<label>>::iterator crfIter =
                cellRemovedFaces.find(fp.masterCell());

            if (crfIter == cellRemovedFaces.end())
            {
                cellRemovedFaces.insert
                (
                    fp.masterCell(),
                    SLList<label>(fp.masterFace())
                );
            }
            else
            {
                crfIter().append(fp.masterFace());
            }

            Map<SLList<face>>::iterator cafIter =
                cellAddedFaces.find(fp.masterCell());
            if (cafIter == cellAddedFaces.end())
            {
                cellAddedFaces.insert
                (
                    fp.masterCell(),
                    SLList<face>(slaveFace.reverseFace())
                );
            }
            else
            {
                cafIter().append(slaveFace.reverseFace());
            }
        }
        else
        {
            // Create cut faces, which replace both master and slave faces

            // Store newly created points
            SLList<point> coupleFacePoints;

            // Master data
            edgeList masterEdges = masterFace.edges();
            List<SLList<label>> masterEdgePoints(masterEdges.size());

            // Slave data
            edgeList slaveEdges = slaveFace.edges();
            List<SLList<label>> slaveEdgePoints(slaveEdges.size());

            // Find common plane
            vector n = masterFace.normal(points_);
            n /= mag(n) + VSMALL;

            // Loop through all edges of the master face. For every edge,
            // intersect it with all edges of the cutting face.
            forAll(masterEdges, masterEdgeI)
            {
                const edge& curMasterEdge = masterEdges[masterEdgeI];

                point P = points_[curMasterEdge.start()];

                // get d and return it into plane
                vector d = curMasterEdge.vec(points_);
                d -= n*(n & d);

                #ifdef DEBUG_COUPLE_INTERSECTION
                Info<< "curMasterEdge: " << curMasterEdge << endl
                    << "P: " << P << endl << "d: " << d << endl;
                #endif

                // go through all slave edges and try to get an intersection.
                // The point is created along the original master edge rather
                // than its corrected direction.
                forAll(slaveEdges, slaveEdgeI)
                {
                    const edge& curSlaveEdge = slaveEdges[slaveEdgeI];

                    point S = points_[curSlaveEdge.start()];

                    // get e and return it into plane
                    vector e = curSlaveEdge.vec(points_);
                    e -= n*(n & e);
                    scalar det = -(e & (n ^ d));

                    #ifdef DEBUG_COUPLE_INTERSECTION
                    Info<< "curSlaveEdge: " << curSlaveEdge << endl
                        << "S: " << S << endl
                        << "e: " << e << endl;
                    #endif

                    if (mag(det) > SMALL)
                    {
                        // non-singular matrix. Look for intersection
                        scalar beta = ((S - P) & (n ^ d))/det;

                        #ifdef DEBUG_COUPLE_INTERSECTION
                        Info<< " beta: " << beta << endl;
                        #endif

                        if (beta > -smallMergeTol_ && beta < 1 + smallMergeTol_)
                        {
                            // slave intersection OK. Try master intersection
                            scalar alpha =
                                (((S - P) & d) + beta*(d & e))/magSqr(d);

                            #ifdef DEBUG_COUPLE_INTERSECTION
                            Info<< " alpha: " << alpha << endl;
                            #endif

                            if
                            (
                                alpha > -smallMergeTol_
                             && alpha < 1 + smallMergeTol_
                            )
                            {
                                // intersection of non-parallel edges
                                #ifdef DEBUG_COUPLE_INTERSECTION
                                Info<< "intersection of non-parallel edges"
                                    << endl;
                                #endif


                                // check for insertion of start-end
                                // points in the middle of the other
                                // edge
                                if (alpha < smallMergeTol_)
                                {
                                    // inserting the start of master edge
                                    if
                                    (
                                        beta > smallMergeTol_
                                     && beta < 1 - smallMergeTol_
                                    )
                                    {
                                        slaveEdgePoints[slaveEdgeI].append
                                        (
                                            curMasterEdge.start()
                                        );
                                    }
                                }
                                else if (alpha > 1 - smallMergeTol_)
                                {
                                    // inserting the end of master edge
                                    if
                                    (
                                        beta > smallMergeTol_
                                     && beta < 1 - smallMergeTol_
                                    )
                                    {
                                        slaveEdgePoints[slaveEdgeI].append
                                        (
                                            curMasterEdge.end()
                                        );
                                    }
                                }
                                else if (beta < smallMergeTol_)
                                {
                                    // inserting the start of the slave edge
                                    if
                                    (
                                        alpha > smallMergeTol_
                                     && alpha < 1 - smallMergeTol_
                                    )
                                    {
                                        masterEdgePoints[masterEdgeI].append
                                        (
                                            curSlaveEdge.start()
                                        );
                                    }
                                }
                                else if (beta > 1 - smallMergeTol_)
                                {
                                    // inserting the start of the slave edge
                                    if
                                    (
                                        alpha > smallMergeTol_
                                     && alpha < 1 - smallMergeTol_
                                    )
                                    {
                                        masterEdgePoints[masterEdgeI].append
                                        (
                                            curSlaveEdge.end()
                                        );
                                    }
                                }
                                else
                                {
                                    masterEdgePoints[masterEdgeI].append
                                    (
                                        nLivePoints + coupleFacePoints.size()
                                    );

                                    slaveEdgePoints[slaveEdgeI].append
                                    (
                                        nLivePoints + coupleFacePoints.size()
                                    );

                                    #ifdef DEBUG_COUPLE_INTERSECTION
                                    Info<< "regular intersection. "
                                        << "Adding point: "
                                        << coupleFacePoints.size()
                                        << " which is "
                                        << P + alpha*curMasterEdge.vec(points_)
                                        << endl;
                                    #endif

                                    // A new point is created. Warning:
                                    // using original edge for accuracy.
                                    //
                                    coupleFacePoints.append
                                        (P + alpha*curMasterEdge.vec(points_));
                                }
                            }
                        }
                    }
                    else
                    {
                        // Add special cases, for intersection of two
                        // parallel line Warning. Here, typically, no new
                        // points will be created. Either one or two of
                        // the slave edge points need to be added to the
                        // master edge and vice versa. The problem is that
                        // no symmetry exists, i.e. both operations needs
                        // to be done separately for both master and slave
                        // side.

                        // Master side
                        // check if the first or second point of slave edge is
                        // on the master edge
                        vector ps = S - P;

                        bool colinear = false;

                        if (mag(ps) < SMALL)
                        {
                            // colinear because P and S are the same point
                            colinear = true;
                        }
                        else if
                        (
                            (ps & d)/(mag(ps)*mag(d)) > 1.0 - smallMergeTol_
                        )
                        {
                            // colinear because ps and d are parallel
                            colinear = true;
                        }

                        if (colinear)
                        {
                            scalar alpha1 = (ps & d)/magSqr(d);

                            if
                            (
                                alpha1 > -smallMergeTol_
                             && alpha1 < 1 + smallMergeTol_
                            )
                            {
                                #ifdef DEBUG_COUPLE_INTERSECTION
                                Info<< "adding irregular master "
                                    << "intersection1: "
                                    << points_[slaveEdges[slaveEdgeI].start()]
                                    << endl;
                                #endif

                                masterEdgePoints[masterEdgeI].append
                                (
                                    slaveEdges[slaveEdgeI].start()
                                );
                            }

                             scalar alpha2 = ((ps + e) & d)/magSqr(d);

                            if
                            (
                                alpha2 > -smallMergeTol_
                             && alpha2 < 1 + smallMergeTol_
                            )
                            {
                                #ifdef DEBUG_COUPLE_INTERSECTION
                                Info<< "adding irregular master "
                                    << "intersection2: "
                                    << points_[slaveEdges[slaveEdgeI].end()]
                                    << endl;
                                #endif

                                masterEdgePoints[masterEdgeI].append
                                (
                                    slaveEdges[slaveEdgeI].end()
                                );
                            }

                            // Slave side
                            // check if the first or second point of
                            // master edge is on the slave edge

                            vector sp = P - S;

                            scalar beta1 = (sp & e)/magSqr(e);

                            #ifdef DEBUG_COUPLE_INTERSECTION
                            Info<< "P: " << P << " S: " << S << " d: " << d
                                << " e: " << e << " sp: " << sp
                                << " beta1: " << beta1 << endl;
                            #endif

                            if
                            (
                                beta1 > -smallMergeTol_
                             && beta1 < 1 + smallMergeTol_
                            )
                            {
                                #ifdef DEBUG_COUPLE_INTERSECTION
                                Info<< "adding irregular slave "
                                    << "intersection1: "
                                    << points_[masterEdges[masterEdgeI].start()]
                                    << endl;
                                #endif

                                slaveEdgePoints[slaveEdgeI].append
                                (
                                    masterEdges[masterEdgeI].start()
                                );
                            }

                            scalar beta2 = ((sp + d) & e)/magSqr(e);

                            if
                            (
                                beta2 > -smallMergeTol_
                             && beta2 < 1 + smallMergeTol_
                            )
                            {
                                #ifdef DEBUG_COUPLE_INTERSECTION
                                Info<< "adding irregular slave "
                                    << "intersection2: "
                                    << points_[masterEdges[masterEdgeI].end()]
                                    << endl;
                                #endif

                                slaveEdgePoints[slaveEdgeI].append
                                (
                                    masterEdges[masterEdgeI].end()
                                );
                            }
                        } // end of colinear
                    } // end of singular intersection
                } // end of slave edges
            } // end of master edges

            #ifdef DEBUG_COUPLE_INTERSECTION
            Info<< "additional slave edge points: " << endl;
            forAll(slaveEdgePoints, edgeI)
            {
                Info<< "edge: " << edgeI << ": " << slaveEdgePoints[edgeI]
                    << endl;
            }
            #endif

            // Add new points
            if (nLivePoints + coupleFacePoints.size() >= points_.size())
            {
                // increase the size of the points list
                Info<< "Resizing points list" << endl;
                points_.setSize(points_.size() + couples_.size());
            }

            for
            (
                SLList<point>::iterator coupleFacePointsIter =
                    coupleFacePoints.begin();
                coupleFacePointsIter != coupleFacePoints.end();
                ++coupleFacePointsIter
            )
            {
                points_[nLivePoints] = coupleFacePointsIter();
                nLivePoints++;
            }

            // edge intersection finished

            // Creating new master side

            // count the number of additional points for face
            label nAdditionalMasterPoints = 0;

            forAll(masterEdgePoints, edgeI)
            {
                nAdditionalMasterPoints += masterEdgePoints[edgeI].size();
            }

            face tmpMasterFace
            (
                masterFace.size()
              + nAdditionalMasterPoints
            );
            label nTmpMasterLabels = 0;

            #ifdef DEBUG_COUPLE_INTERSECTION
            Info<< "masterFace: " << masterFace << endl
                << "nAdditionalMasterPoints: " << nAdditionalMasterPoints
                << endl;
            #endif

            forAll(masterEdges, masterEdgeI)
            {
                // Insert the starting point of the edge
                tmpMasterFace[nTmpMasterLabels] =
                    masterEdges[masterEdgeI].start();
                nTmpMasterLabels++;

                // get reference to added points of current edge
                const SLList<label>& curMEdgePoints =
                    masterEdgePoints[masterEdgeI];

                // create a markup list of points that have been used
                boolList usedMasterPoint(curMEdgePoints.size(), false);

                vector edgeVector = masterEdges[masterEdgeI].vec(points_);

                #ifdef DEBUG_FACE_ORDERING
                Info<< "edgeVector: " << edgeVector << endl
                    << "curMEdgePoints.size(): " << curMEdgePoints.size()
                    << endl;
                #endif

                // renormalise
                edgeVector /= magSqr(edgeVector);

                point edgeStartPoint =
                    points_[masterEdges[masterEdgeI].start()];

                // loop until the next label to add is -1
                for (;;)
                {
                    label nextPointLabel = -1;
                    label usedI = -1;
                    scalar minAlpha = GREAT;

                    label i = 0;

                    for
                    (
                        SLList<label>::const_iterator curMEdgePointsIter =
                            curMEdgePoints.begin();
                        curMEdgePointsIter != curMEdgePoints.end();
                        ++curMEdgePointsIter
                    )
                    {
                        if (!usedMasterPoint[i])
                        {
                            scalar alpha =
                                edgeVector
                              & (
                                    points_[curMEdgePointsIter()]
                                  - edgeStartPoint
                                );

                            #ifdef DEBUG_FACE_ORDERING
                            Info<< " edgeStartPoint: " << edgeStartPoint
                                << " edgeEndPoint: "
                                << points_[masterEdges[masterEdgeI].end()]
                                << " other point: "
                                << points_[curMEdgePointsIter()]
                                << " alpha: " << alpha << endl;
                            #endif

                            if (alpha < minAlpha)
                            {
                                minAlpha = alpha;
                                usedI = i;
                                nextPointLabel = curMEdgePointsIter();
                            }
                        }

                        #ifdef DEBUG_FACE_ORDERING
                        Info<< "nextPointLabel: " << nextPointLabel << endl;
                        #endif

                        i++;
                    }

                    if (nextPointLabel > -1)
                    {
                        #ifdef DEBUG_FACE_ORDERING
                        Info<< "added nextPointLabel: " << nextPointLabel
                            << " nTmpMasterLabels: " << nTmpMasterLabels
                            << " to place " << nTmpMasterLabels << endl;
                        #endif

                        usedMasterPoint[usedI] = true;
                        // add the next point
                        tmpMasterFace[nTmpMasterLabels] =
                            nextPointLabel;
                        nTmpMasterLabels++;
                    }
                    else
                    {
                        break;
                    }
                }
            }

            // reset the size of master
            tmpMasterFace.setSize(nTmpMasterLabels);

            #ifdef DEBUG_FACE_ORDERING
            Info<< "tmpMasterFace: " << tmpMasterFace << endl;
            #endif

            // Eliminate all zero-length edges
            face newMasterFace(labelList(tmpMasterFace.size(), labelMax));

            // insert first point by hand. Careful: the first one is
            // used for comparison to allow the edge collapse across
            // point zero.
            //
            newMasterFace[0] = tmpMasterFace[0];
            label nMaster = 0;

            edgeList mstEdgesToCollapse = tmpMasterFace.edges();

            scalar masterTol =
                cpMergePointTol_*boundBox(tmpMasterFace.points(points_)).mag();

            forAll(mstEdgesToCollapse, edgeI)
            {
                #ifdef DEBUG_FACE_ORDERING
                Info<< "edgeI: " << edgeI << " curEdge: "
                    << mstEdgesToCollapse[edgeI] << endl
                    << "master edge " << edgeI << ", "
                    << mstEdgesToCollapse[edgeI].mag(points_) << endl;
                #endif

                // Edge merge tolerance = masterTol
                if (mstEdgesToCollapse[edgeI].mag(points_) < masterTol)
                {
                    newMasterFace[nMaster] =
                        min
                        (
                            newMasterFace[nMaster],
                            mstEdgesToCollapse[edgeI].end()
                        );

                    #ifdef DEBUG_FACE_ORDERING
                    Info<< "Collapsed: nMaster: " << nMaster
                        << " label: " << newMasterFace[nMaster] << endl;
                    #endif

                }
                else
                {
                    nMaster++;

                    if (edgeI < mstEdgesToCollapse.size() - 1)
                    {
                        // last edge does not add the point
                    #ifdef DEBUG_FACE_ORDERING
                        Info<< "Added: nMaster: " << nMaster
                            << " label: " << mstEdgesToCollapse[edgeI].end()
                            << endl;
                    #endif

                        newMasterFace[nMaster] =
                            mstEdgesToCollapse[edgeI].end();
                    }
                }
            }

            newMasterFace.setSize(nMaster);

            #ifdef DEBUG_COUPLE
            Info<< "newMasterFace: " << newMasterFace << endl
                << "points: " << newMasterFace.points(points_) << endl;
            #endif

            // Creating new slave side

            // count the number of additional points for face
            label nAdditionalSlavePoints = 0;

            forAll(slaveEdgePoints, edgeI)
            {
                nAdditionalSlavePoints += slaveEdgePoints[edgeI].size();
            }

            face tmpSlaveFace
            (
                slaveFace.size()
              + nAdditionalSlavePoints
            );
            label nTmpSlaveLabels = 0;

            #ifdef DEBUG_COUPLE_INTERSECTION
            Info<< "slaveFace: " << slaveFace << endl
                << "nAdditionalSlavePoints: " << nAdditionalSlavePoints << endl;
            #endif

            forAll(slaveEdges, slaveEdgeI)
            {
                // Insert the starting point of the edge
                tmpSlaveFace[nTmpSlaveLabels] =
                    slaveEdges[slaveEdgeI].start();
                nTmpSlaveLabels++;

                // get reference to added points of current edge
                const SLList<label>& curSEdgePoints =
                    slaveEdgePoints[slaveEdgeI];

                // create a markup list of points that have been used
                boolList usedSlavePoint(curSEdgePoints.size(), false);

                vector edgeVector = slaveEdges[slaveEdgeI].vec(points_);

                #ifdef DEBUG_FACE_ORDERING
                Info<< "curSEdgePoints.size(): "
                    << curSEdgePoints.size() << endl
                    << "edgeVector: " << edgeVector << endl;
                #endif

                // renormalise
                edgeVector /= magSqr(edgeVector);

                point edgeStartPoint =
                    points_[slaveEdges[slaveEdgeI].start()];

                // loop until the next label to add is -1
                for (;;)
                {
                    label nextPointLabel = -1;
                    label usedI = -1;
                    scalar minAlpha = GREAT;

                    label i = 0;

                    for
                    (
                        SLList<label>::const_iterator curSEdgePointsIter =
                            curSEdgePoints.begin();
                        curSEdgePointsIter != curSEdgePoints.end();
                        ++curSEdgePointsIter
                    )
                    {
                        if (!usedSlavePoint[i])
                        {
                            scalar alpha =
                                edgeVector
                              & (
                                    points_[curSEdgePointsIter()]
                                  - edgeStartPoint
                                );

                            #ifdef DEBUG_FACE_ORDERING
                            Info<< " edgeStartPoint: " << edgeStartPoint
                                << " edgeEndPoint: "
                                << points_[slaveEdges[slaveEdgeI].end()]
                                << " other point: "
                                << points_[curSEdgePointsIter()]
                                << " alpha: " << alpha << endl;
                            #endif

                            if (alpha < minAlpha)
                            {
                                minAlpha = alpha;
                                usedI = i;
                                nextPointLabel = curSEdgePointsIter();
                            }
                        }

                        #ifdef DEBUG_FACE_ORDERING
                        Info<< "nextPointLabel: " << nextPointLabel << endl;
                        #endif

                        i++;
                    }

                    if (nextPointLabel > -1)
                    {
                        #ifdef DEBUG_FACE_ORDERING
                        Info<< "added nextPointLabel: " << nextPointLabel
                            << " nTmpSlaveLabels: " << nTmpSlaveLabels
                            << " to place " << nTmpSlaveLabels << endl;
                        #endif

                        usedSlavePoint[usedI] = true;
                        // add the next point
                        tmpSlaveFace[nTmpSlaveLabels] =
                            nextPointLabel;
                        nTmpSlaveLabels++;
                    }
                    else
                    {
                        break;
                    }
                }
            }

            // reset the size of slave
            tmpSlaveFace.setSize(nTmpSlaveLabels);

            #ifdef DEBUG_FACE_ORDERING
            Info<< "tmpSlaveFace: " << tmpSlaveFace << endl;
            #endif

            // Eliminate all zero-length edges
            face newSlaveFace(labelList(tmpSlaveFace.size(), labelMax));

            // insert first point by hand. Careful: the first one is
            // used for comparison to allow the edge collapse across
            // point zero.
            //
            newSlaveFace[0] = tmpSlaveFace[0];
            label nSlave = 0;

            edgeList slvEdgesToCollapse = tmpSlaveFace.edges();

            scalar slaveTol =
                cpMergePointTol_*boundBox(tmpSlaveFace.points(points_)).mag();

            forAll(slvEdgesToCollapse, edgeI)
            {
                #ifdef DEBUG_FACE_ORDERING
                Info<< "slave edge length: " << edgeI << ", "
                    << slvEdgesToCollapse[edgeI].mag(points_)<< endl;
                #endif

                 // edge merge tolerance = slaveTol
                if (slvEdgesToCollapse[edgeI].mag(points_) < slaveTol)
                {
                    newSlaveFace[nSlave] =
                        min
                        (
                            newSlaveFace[nSlave],
                            slvEdgesToCollapse[edgeI].end()
                        );
                }
                else
                {
                    nSlave++;
                    if (edgeI < slvEdgesToCollapse.size() - 1)
                    {
                        // last edge does not add the point
                        newSlaveFace[nSlave] = slvEdgesToCollapse[edgeI].end();
                    }
                }
            }

            newSlaveFace.setSize(nSlave);

            #ifdef DEBUG_COUPLE
            Info<< "newSlaveFace: " << newSlaveFace << endl
                << "points: " << newSlaveFace.points(points_) << endl << endl;
            #endif

            // Create the intersection face

            // Algorithm:
            // Loop through
            // points of the master and try to find one which falls
            // within the slave.  If not found, look through all
            // edges of the slave and find one which falls within the
            // master.  This point will be the starting location for
            // the cut face.

            edgeList newMasterEdges = newMasterFace.edges();
            edgeList newSlaveEdges = newSlaveFace.edges();

            #ifdef DEBUG_RIGHT_HAND_WALK
            Info<< "newMasterEdges: " << newMasterEdges << endl
                << "newSlaveEdges: " << newSlaveEdges << endl;
            #endif

            edge startEdge(-1, -1);

            // Remember where the start edge was found:
            // 0 for not found
            // 1 for master
            // 2 for slave
            label startEdgeFound = 0;

            vector masterProjDir = -newMasterFace.normal(points_);

            forAll(newSlaveEdges, edgeI)
            {
                // Take the slave edge points and project into the master.
                // In order to create a good intersection, move the
                // point away from the master in the direction of its
                // normal.
                point pointStart = points_[newSlaveEdges[edgeI].start()];

                point pointEnd = points_[newSlaveEdges[edgeI].end()];

                if
                (
                    newMasterFace.ray
                    (
                        pointStart,
                        masterProjDir,
                        points_,
                        intersection::FULL_RAY
                    ).hit()
                 && newMasterFace.ray
                    (
                        pointEnd,
                        masterProjDir,
                        points_,
                        intersection::FULL_RAY
                    ).hit()
                )
                {
                    startEdge = newSlaveEdges[edgeI];
                    startEdgeFound = 2;

                    #ifdef DEBUG_RIGHT_HAND_WALK
                    Info<< "slave edge found" << endl;
                    #endif

                    break;
                }
            }

            if (startEdgeFound == 0)
            {
                vector slaveProjDir = -newSlaveFace.normal(points_);

                forAll(newMasterEdges, edgeI)
                {
                    // Take the edge master points and project into the slave.
                    // In order to create a good intersection, move the
                    // point away from the slave in the direction of its
                    // normal.
                    point pointStart = points_[newMasterEdges[edgeI].start()];

                    point pointEnd = points_[newMasterEdges[edgeI].end()];

                    if
                    (
                        newSlaveFace.ray
                        (
                            pointStart,
                            slaveProjDir,
                            points_,
                        intersection::FULL_RAY
                        ).hit()
                     && newSlaveFace.ray
                        (
                            pointEnd,
                            slaveProjDir,
                            points_,
                            intersection::FULL_RAY
                        ).hit()
                    )
                    {
                        startEdge = newMasterEdges[edgeI];
                        startEdgeFound = 1;

                        #ifdef DEBUG_RIGHT_HAND_WALK
                        Info<< "master edge found" << endl;
                        #endif

                        break;
                    }
                }
            }

            // create the intersected face using right-hand walk rule
            face intersectedFace
            (
                labelList(newMasterFace.size() + newSlaveFace.size(), -1)
            );

            if (startEdgeFound > 0)
            {
                #ifdef DEBUG_RIGHT_HAND_WALK
                Info<< "start edge: " << startEdge << endl;
                #endif

                // Loop through both faces and add all edges
                // containing the current point and add them to the
                // list of edges to consider.  Make sure all edges are
                // added such that the current point is their start.
                // Loop through all edges to consider and find the one
                // which produces the buggest right-hand-turn.  This
                // is the next edge to be added to the face.  If its
                // end is the same as the starting point, the face is
                // complete; resize it to the number of active points
                // and exit.

                vector planeNormal = newMasterFace.normal(points_);
                planeNormal /= mag(planeNormal) + VSMALL;

                #ifdef DEBUG_RIGHT_HAND_WALK
                Info<< "planeNormal: " << planeNormal << endl;
                #endif

                // Do a check to control the right-hand turn.  This is
                // based on the triple product of the edge start
                // vector to face centre, the edge vector and the
                // plane normal.  If the triple product is negative,
                // the edge needs to be reversed to allow the
                // right-hand-turn rule to work.

                vector faceCentre;

                if (startEdgeFound == 1)
                {
                    faceCentre = newMasterFace.centre(points_);
                }
                else
                {
                    faceCentre = newSlaveFace.centre(points_);
                }

                scalar tripleProduct =
                    (
                        (faceCentre - points_[startEdge.start()])
                      ^ startEdge.vec(points_)
                    ) & planeNormal;

                if (tripleProduct < 0)
                {
                    #ifdef DEBUG_RIGHT_HAND_WALK
                    Info<< "Turning edge for right-hand turn rule" << endl;
                    #endif
                    startEdge.flip();
                }

                // prepare the loop for the right-hand walk
                intersectedFace[0] = startEdge.start();
                intersectedFace[1] = startEdge.end();
                label nIntFacePoints = 2;

                edge curEdge = startEdge;

                bool completedFace = false;

                do
                {
                    SLList<edge> edgesToConsider;

                    // collect master edges
                    forAll(newMasterEdges, edgeI)
                    {
                        const edge& cme = newMasterEdges[edgeI];

                        if (cme != curEdge)
                        {
                            if (cme.start() == curEdge.end())
                            {
                                edgesToConsider.append(cme);
                            }
                            else if (cme.end() == curEdge.end())
                            {
                                edgesToConsider.append(cme.reverseEdge());
                            }
                            // otherwise, it does not have the current point
                        }
                    }

                    // collect slave edges
                    forAll(newSlaveEdges, edgeI)
                    {
                        const edge& cse = newSlaveEdges[edgeI];

                        if (cse != curEdge)
                        {
                            if (cse.start() == curEdge.end())
                            {
                                edgesToConsider.append(cse);
                            }
                            else if (cse.end() == curEdge.end())
                            {
                                edgesToConsider.append(cse.reverseEdge());
                            }
                            // otherwise, it does not have the current point
                        }
                    }

                    #ifdef DEBUG_RIGHT_HAND_WALK
                    Info<< "number of edges to consider: "
                        << edgesToConsider.size() << endl
                        << "edges to consider: " << edgesToConsider << endl;
                    #endif

                    if (edgesToConsider.empty())
                    {
                        FatalErrorInFunction
                            << setprecision(12)
                            << "void starMesh::createCoupleMatches() : "
                            << endl << "error in face intersection: "
                            << "no edges to consider for closing the loop"
                            << coupleI << ". STAR couple ID: "
                            << couples_[coupleI].coupleID() << endl
                            << "Cut Master Face: " << newMasterFace << endl
                            << "points: " << newMasterFace.points(points_)
                            << endl
                            << "Cut Slave Face: " << newSlaveFace << endl
                            << "points: " << newSlaveFace.points(points_)
                            << endl << "intersected face: "
                            << intersectedFace
                            << abort(FatalError);
                    }

                    // vector along the edge
                    vector ahead = curEdge.vec(points_);
                    ahead -= planeNormal*(planeNormal & ahead);
                    ahead /= mag(ahead) + VSMALL;

                    // vector pointing right
                    vector right = ahead ^ planeNormal;
                    right /= mag(right) + VSMALL;

                    // first edge taken for reference
                    edge nextEdge = edgesToConsider.first();
                    vector nextEdgeVec = nextEdge.vec(points_);
                    nextEdgeVec -= planeNormal*(planeNormal & nextEdgeVec);
                    nextEdgeVec /= mag(nextEdgeVec) + VSMALL;

                    scalar rightTurn = nextEdgeVec & right;
                    scalar goStraight = nextEdgeVec & ahead;

                    #ifdef DEBUG_RIGHT_HAND_WALK
                    Info<< "rightTurn: " << rightTurn
                        << " goStraight: " << goStraight << endl;
                    #endif

                    for
                    (
                        SLList<edge>::iterator etcIter =
                            edgesToConsider.begin();
                        etcIter != edgesToConsider.end();
                        ++etcIter
                    )
                    {
                        // right-hand walk rule
                        vector newDir = etcIter().vec(points_);
                        newDir -= planeNormal*(planeNormal & newDir);
                        newDir /= mag(newDir) + VSMALL;

                        scalar curRightTurn = newDir & right;
                        scalar curGoStraight = newDir & ahead;

                        #ifdef DEBUG_RIGHT_HAND_WALK
                        Info<< "curRightTurn: " << curRightTurn
                            << " curGoStraight: " << curGoStraight << endl;
                        #endif

                        if (rightTurn < 0) // old edge turning left
                        {
                            if (curRightTurn < 0) // new edge turning left
                            {
                                // both go left. Grab one with greater ahead
                                if (curGoStraight > goStraight)
                                {
                                    #ifdef DEBUG_RIGHT_HAND_WALK
                                    Info<< "a" << endl;
                                    #endif

                                    // Good edge, turning left less than before
                                    nextEdge = etcIter();
                                    rightTurn = curRightTurn;
                                    goStraight = curGoStraight;
                                }
                            }
                            else // new edge turning right
                            {
                                #ifdef DEBUG_RIGHT_HAND_WALK
                                Info<< "b" << endl;
                                #endif

                                // good edge, turning right
                                nextEdge = etcIter();
                                rightTurn = curRightTurn;
                                goStraight = curGoStraight;
                            }
                        }
                        else // old edge turning right
                        {
                            // new edge turning left rejected
                            if (curRightTurn >= 0) // new edge turning right
                            {
                                // grab one with smaller ahead
                                if (curGoStraight < goStraight)
                                {
                                    #ifdef DEBUG_RIGHT_HAND_WALK
                                    Info<< "c" << endl;
                                    #endif

                                    // good edge, turning right more than before
                                    nextEdge = etcIter();
                                    rightTurn = curRightTurn;
                                    goStraight = curGoStraight;
                                }
                            }
                        }
                    }

                    // check if the loop is completed
                    if (nextEdge.end() == intersectedFace[0])
                    {
                        // loop is completed. No point to add
                        completedFace = true;
                    }
                    else
                    {
                        // Check if there is room for more points
                        if (nIntFacePoints >= intersectedFace.size())
                        {
                            FatalErrorInFunction
                                << setprecision(12)
                                << "void starMesh::createCoupleMatches() : "
                                << endl << "error in intersected face: "
                                << "lost thread for intersection of couple "
                                << coupleI << ". STAR couple ID: "
                                << couples_[coupleI].coupleID() << endl
                                << "Cut Master Face: " << newMasterFace << endl
                                << "points: " << newMasterFace.points(points_)
                                << endl
                                << "Cut Slave Face: " << newSlaveFace << endl
                                << "points: " << newSlaveFace.points(points_)
                                << endl << "intersected face: "
                                << intersectedFace
                                << abort(FatalError);
                        }

                        // insert the point
                        intersectedFace[nIntFacePoints] = nextEdge.end();
                        nIntFacePoints++;

                        // grab the current point and the current edge
                        curEdge = nextEdge;

                        #ifdef DEBUG_RIGHT_HAND_WALK
                        Info<< "inserted point " << nextEdge.end() << endl
                            << "curEdge: " << curEdge << endl;
                        #endif
                    }
                }
                while (!completedFace);

                // resize the face
                intersectedFace.setSize(nIntFacePoints);

                #ifdef DEBUG_COUPLE
                Info<< "intersectedFace: " << intersectedFace << endl;
                #endif

                // check the intersection face for duplicate points
                forAll(intersectedFace, checkI)
                {
                    for
                    (
                        label checkJ = checkI + 1;
                        checkJ < intersectedFace.size();
                        checkJ++
                    )
                    {
                        if (intersectedFace[checkI] == intersectedFace[checkJ])
                        {
                            FatalErrorInFunction
                                << setprecision(12)
                                << "void starMesh::createCoupleMatches() : "
                                << endl << "error in intersected face: "
                                << "duplicate point in intersected face "
                                << "for couple no " << coupleI
                                << ". STAR couple ID: "
                                << couples_[coupleI].coupleID() << endl
                                << "Duplicate point label: "
                                << intersectedFace[checkI] << endl
                                << "Cut Master Face: " << newMasterFace << endl
                                << "points: " << newMasterFace.points(points_)
                                << endl
                                << "Cut Slave Face: " << newSlaveFace << endl
                                << "points: " << newSlaveFace.points(points_)
                                << endl << "intersected face: "
                                << intersectedFace
                                << abort(FatalError);
                        }
                    }
                }
            }
            else
            {
                FatalErrorInFunction
                    << setprecision(12)
                    << "void starMesh::createCoupleMatches() : " << endl
                    << "could not find start edge for intersection of couple "
                    << coupleI << ". STAR couple ID: "
                    << couples_[coupleI].coupleID() << endl
                    << "Cut Master Face: " << newMasterFace << endl
                    << "points: " << newMasterFace.points(points_) << endl
                    << "Cut Slave Face: " << newSlaveFace << endl
                    << "points: " << newSlaveFace.points(points_)
                    << abort(FatalError);
            }

            // Project all points of the intersected face
            // onto the master face to ensure closedness
            vector pointProjectionNormal = -masterFace.normal(points_);

            forAll(intersectedFace, intPointi)
            {
                #ifdef DEBUG_COUPLE_PROJECTION
                Info<< "Proj: old point: "
                    << points_[intersectedFace[intPointi]] << endl;
                #endif

                pointHit projHit =
                    masterFace.ray
                    (
                        points_[intersectedFace[intPointi]],
                        pointProjectionNormal,
                        points_,
                        intersection::FULL_RAY
                    );

                if (projHit.hit())
                {
                    points_[intersectedFace[intPointi]] =
                        projHit.hitPoint();

                    #ifdef DEBUG_COUPLE_PROJECTION
                    Info<< "      new point: "
                        << points_[intersectedFace[intPointi]] << endl;
                    #endif
                }
            }

            // Check the direction of the intersection face
            if
            (
                (
                    masterFace.normal(points_)
                  & intersectedFace.normal(points_)
                ) < VSMALL
            )
            {
                intersectedFace.flip();
            }

            // Add the new face to both master and slave

            // Master face is replaced by a set of slave faces
            Map<SLList<label>>::iterator crfMasterIter =
                cellRemovedFaces.find(fp.masterCell());

            if (crfMasterIter == cellRemovedFaces.end())
            {
                cellRemovedFaces.insert
                (
                    fp.masterCell(),
                    SLList<label>(fp.masterFace())
                );
            }
            else
            {
                crfMasterIter().append(fp.masterFace());
            }

            Map<SLList<label>>::iterator crfSlaveIter =
                cellRemovedFaces.find(fp.slaveCell());

            if (crfSlaveIter == cellRemovedFaces.end())
            {
                cellRemovedFaces.insert
                (
                    fp.slaveCell(),
                    SLList<label>(fp.slaveFace())
                );
            }
            else
            {
                crfSlaveIter().append(fp.slaveFace());
            }

            Map<SLList<face>>::iterator cafMasterIter =
                cellAddedFaces.find(fp.masterCell());
            if (cafMasterIter == cellAddedFaces.end())
            {
                cellAddedFaces.insert
                (
                    fp.masterCell(),
                    SLList<face>(intersectedFace)
                );
            }
            else
            {
                cafMasterIter().append(intersectedFace);
            }

            Map<SLList<face>>::iterator cafSlaveIter =
                cellAddedFaces.find(fp.slaveCell());
            if (cafSlaveIter == cellAddedFaces.end())
            {
                cellAddedFaces.insert
                (
                    fp.slaveCell(),
                    SLList<face>(intersectedFace.reverseFace())
                );
            }
            else
            {
                cafSlaveIter().append(intersectedFace.reverseFace());
            }
        } // end of arbitrary match
    }

    if (couples_.size())
    {
        // Loop through all cells and reset faces for removal to zero size
        const labelList crfToc = cellRemovedFaces.toc();

        forAll(crfToc, celli)
        {
            const label curCell = crfToc[celli];

            const SLList<label>& curRemovedFaces = cellRemovedFaces[curCell];

            for
            (
                SLList<label>::const_iterator curRemovedFacesIter =
                    curRemovedFaces.begin();
                curRemovedFacesIter != curRemovedFaces.end();
                ++curRemovedFacesIter
            )
            {
                cellFaces_[curCell][curRemovedFacesIter()].setSize(0);
            }

            if (curRemovedFaces.size())
            {
                // reset the shape pointer to unknown
                cellShapes_[curCell] = cellShape(*unknownPtr_, labelList(0));
            }
        }

        const labelList cafToc = cellAddedFaces.toc();

        // Insert the new faces into the list
        forAll(cafToc, celli)
        {
            const label curCell = cafToc[celli];

            const SLList<face>& curAddedFaces = cellAddedFaces[curCell];

            faceList oldFaces = cellFaces_[curCell];

            faceList& newFaces = cellFaces_[curCell];

            newFaces.setSize(oldFaces.size() + curAddedFaces.size());
            label nNewFaces = 0;

            // copy original faces that have not been removed
            forAll(oldFaces, facei)
            {
                if (oldFaces[facei].size())
                {
                    newFaces[nNewFaces] = oldFaces[facei];
                    nNewFaces++;
                }
            }

            // add new faces
            for
            (
                SLList<face>::const_iterator curAddedFacesIter =
                    curAddedFaces.begin();
                curAddedFacesIter != curAddedFaces.end();
                ++curAddedFacesIter
            )
            {
                newFaces[nNewFaces] = curAddedFacesIter();
                nNewFaces++;
            }

            // reset the size of the face list
            newFaces.setSize(nNewFaces);

            if (curAddedFaces.size())
            {
                // reset the shape pointer to unknown
                cellShapes_[curCell] = cellShape(*unknownPtr_, labelList(0));
            }
        }

        // Resize the point list to the number of created points
        points_.setSize(nLivePoints);

        // Finished
        Info<< "Finished doing couples" << endl;
    }
}


// ************************************************************************* //