File: meshDualiser.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 (1439 lines) | stat: -rw-r--r-- 41,679 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
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  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/>.

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

#include "meshDualiser.H"
#include "meshTools.H"
#include "polyMesh.H"
#include "polyTopoChange.H"
#include "mapPolyMesh.H"
#include "edgeFaceCirculator.H"
#include "mergePoints.H"
#include "OFstream.H"

// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //

namespace Foam
{
    defineTypeNameAndDebug(meshDualiser, 0);
}


// * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //

void Foam::meshDualiser::checkPolyTopoChange(const polyTopoChange& meshMod)
{
    // Assume no removed points
    pointField points(meshMod.points().size());
    forAll(meshMod.points(), i)
    {
        points[i] = meshMod.points()[i];
    }

    labelList oldToNew;
    label nUnique = mergePoints
    (
        points,
        1e-6,
        false,
        oldToNew
    );

    if (nUnique < points.size())
    {
        labelListList newToOld(invertOneToMany(nUnique, oldToNew));

        forAll(newToOld, newI)
        {
            if (newToOld[newI].size() != 1)
            {
                FatalErrorInFunction
                    << "duplicate verts:" << newToOld[newI]
                    << " coords:"
                    << UIndirectList<point>(points, newToOld[newI])()
                    << abort(FatalError);
            }
        }
    }
}


// Dump state so far.
void Foam::meshDualiser::dumpPolyTopoChange
(
    const polyTopoChange& meshMod,
    const fileName& prefix
)
{
    OFstream str1(prefix + "Faces.obj");
    OFstream str2(prefix + "Edges.obj");

    Info<< "Dumping current polyTopoChange. Faces to " << str1.name()
        << " , points and edges to " << str2.name() << endl;

    const DynamicList<point>& points = meshMod.points();

    forAll(points, pointi)
    {
        meshTools::writeOBJ(str1, points[pointi]);
        meshTools::writeOBJ(str2, points[pointi]);
    }

    const DynamicList<face>& faces = meshMod.faces();

    forAll(faces, facei)
    {
        const face& f = faces[facei];

        str1<< 'f';
        forAll(f, fp)
        {
            str1<< ' ' << f[fp]+1;
        }
        str1<< nl;

        str2<< 'l';
        forAll(f, fp)
        {
            str2<< ' ' << f[fp]+1;
        }
        str2<< ' ' << f[0]+1 << nl;
    }
}


Foam::label Foam::meshDualiser::findDualCell
(
    const label celli,
    const label pointi
) const
{
    const labelList& dualCells = pointToDualCells_[pointi];

    if (dualCells.size() == 1)
    {
        return dualCells[0];
    }
    else
    {
        label index = findIndex(mesh_.pointCells()[pointi], celli);

        return dualCells[index];
    }
}


void Foam::meshDualiser::generateDualBoundaryEdges
(
    const PackedBoolList& isBoundaryEdge,
    const label pointi,
    polyTopoChange& meshMod
)
{
    const labelList& pEdges = mesh_.pointEdges()[pointi];

    forAll(pEdges, pEdgeI)
    {
        label edgeI = pEdges[pEdgeI];

        if (edgeToDualPoint_[edgeI] == -1 && isBoundaryEdge.get(edgeI) == 1)
        {
            const edge& e = mesh_.edges()[edgeI];

            edgeToDualPoint_[edgeI] = meshMod.addPoint
            (
                e.centre(mesh_.points()),
                pointi, // masterPoint
                -1,     // zoneID
                true    // inCell
            );
        }
    }
}


// Return true if point on face has same dual cells on both owner and neighbour
// sides.
bool Foam::meshDualiser::sameDualCell
(
    const label facei,
    const label pointi
) const
{
    if (!mesh_.isInternalFace(facei))
    {
        FatalErrorInFunction
            << "face:" << facei << " is not internal face."
            << abort(FatalError);
    }

    label own = mesh_.faceOwner()[facei];
    label nei = mesh_.faceNeighbour()[facei];

    return findDualCell(own, pointi) == findDualCell(nei, pointi);
}


Foam::label Foam::meshDualiser::addInternalFace
(
    const label masterPointi,
    const label masterEdgeI,
    const label masterFacei,

    const bool edgeOrder,
    const label dualCell0,
    const label dualCell1,
    const DynamicList<label>& verts,
    polyTopoChange& meshMod
) const
{
    face newFace(verts);

    if (edgeOrder != (dualCell0 < dualCell1))
    {
        reverse(newFace);
    }

    if (debug)
    {
        pointField facePoints(meshMod.points(), newFace);

        labelList oldToNew;
        label nUnique = mergePoints
        (
            facePoints,
            1e-6,
            false,
            oldToNew
        );

        if (nUnique < facePoints.size())
        {
            FatalErrorInFunction
                << "verts:" << verts << " newFace:" << newFace
                << " face points:" << facePoints
                << abort(FatalError);
        }
    }


    label zoneID = -1;
    bool zoneFlip = false;
    if (masterFacei != -1)
    {
        zoneID = mesh_.faceZones().whichZone(masterFacei);

        if (zoneID != -1)
        {
            const faceZone& fZone = mesh_.faceZones()[zoneID];

            zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
        }
    }

    label dualFacei;

    if (dualCell0 < dualCell1)
    {
        dualFacei = meshMod.addFace
        (
            newFace,
            dualCell0,      // own
            dualCell1,      // nei
            masterPointi,   // masterPointID
            masterEdgeI,    // masterEdgeID
            masterFacei,    // masterFaceID
            false,          // flipFaceFlux
            -1,             // patchID
            zoneID,         // zoneID
            zoneFlip        // zoneFlip
        );

        //pointField dualPoints(meshMod.points());
        //vector n(newFace.normal(dualPoints));
        //n /= mag(n);
        //Pout<< "Generated internal dualFace:" << dualFacei
        //    << " verts:" << newFace
        //    << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
        //    << " n:" << n
        //    << " between dualowner:" << dualCell0
        //    << " dualneigbour:" << dualCell1
        //    << endl;
    }
    else
    {
        dualFacei = meshMod.addFace
        (
            newFace,
            dualCell1,      // own
            dualCell0,      // nei
            masterPointi,   // masterPointID
            masterEdgeI,    // masterEdgeID
            masterFacei,    // masterFaceID
            false,          // flipFaceFlux
            -1,             // patchID
            zoneID,         // zoneID
            zoneFlip        // zoneFlip
        );

        //pointField dualPoints(meshMod.points());
        //vector n(newFace.normal(dualPoints));
        //n /= mag(n);
        //Pout<< "Generated internal dualFace:" << dualFacei
        //    << " verts:" << newFace
        //    << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
        //    << " n:" << n
        //    << " between dualowner:" << dualCell1
        //    << " dualneigbour:" << dualCell0
        //    << endl;
    }
    return dualFacei;
}


Foam::label Foam::meshDualiser::addBoundaryFace
(
    const label masterPointi,
    const label masterEdgeI,
    const label masterFacei,

    const label dualCelli,
    const label patchi,
    const DynamicList<label>& verts,
    polyTopoChange& meshMod
) const
{
    face newFace(verts);

    label zoneID = -1;
    bool zoneFlip = false;
    if (masterFacei != -1)
    {
        zoneID = mesh_.faceZones().whichZone(masterFacei);

        if (zoneID != -1)
        {
            const faceZone& fZone = mesh_.faceZones()[zoneID];

            zoneFlip = fZone.flipMap()[fZone.whichFace(masterFacei)];
        }
    }

    label dualFacei = meshMod.addFace
    (
        newFace,
        dualCelli,      // own
        -1,             // nei
        masterPointi,   // masterPointID
        masterEdgeI,    // masterEdgeID
        masterFacei,    // masterFaceID
        false,          // flipFaceFlux
        patchi,         // patchID
        zoneID,         // zoneID
        zoneFlip        // zoneFlip
    );

    //pointField dualPoints(meshMod.points());
    //vector n(newFace.normal(dualPoints));
    //n /= mag(n);
    //Pout<< "Generated boundary dualFace:" << dualFacei
    //    << " verts:" << newFace
    //    << " points:" << UIndirectList<point>(meshMod.points(), newFace)()
    //    << " n:" << n
    //    << " on dualowner:" << dualCelli
    //    << endl;
    return dualFacei;
}


// Walks around edgeI.
// splitFace=true : creates multiple faces
// splitFace=false: creates single face if same dual cells on both sides,
//                  multiple faces otherwise.
void Foam::meshDualiser::createFacesAroundEdge
(
    const bool splitFace,
    const PackedBoolList& isBoundaryEdge,
    const label edgeI,
    const label startFacei,
    polyTopoChange& meshMod,
    boolList& doneEFaces
) const
{
    const edge& e = mesh_.edges()[edgeI];
    const labelList& eFaces = mesh_.edgeFaces()[edgeI];

    label fp = edgeFaceCirculator::getMinIndex
    (
        mesh_.faces()[startFacei],
        e[0],
        e[1]
    );

    edgeFaceCirculator ie
    (
        mesh_,
        startFacei, // face
        true,       // ownerSide
        fp,         // fp
        isBoundaryEdge.get(edgeI) == 1  // isBoundaryEdge
    );
    ie.setCanonical();

    bool edgeOrder = ie.sameOrder(e[0], e[1]);
    label startFaceLabel = ie.faceLabel();

    //Pout<< "At edge:" << edgeI << " verts:" << e
    //    << " points:" << mesh_.points()[e[0]] << mesh_.points()[e[1]]
    //    << " started walking at face:" << ie.faceLabel()
    //    << " verts:" << mesh_.faces()[ie.faceLabel()]
    //    << " edgeOrder:" << edgeOrder
    //    << " in direction of cell:" << ie.cellLabel()
    //    << endl;

    // Walk and collect face.
    DynamicList<label> verts(100);

    if (edgeToDualPoint_[edgeI] != -1)
    {
        verts.append(edgeToDualPoint_[edgeI]);
    }
    if (faceToDualPoint_[ie.faceLabel()] != -1)
    {
        doneEFaces[findIndex(eFaces, ie.faceLabel())] = true;
        verts.append(faceToDualPoint_[ie.faceLabel()]);
    }
    if (cellToDualPoint_[ie.cellLabel()] != -1)
    {
        verts.append(cellToDualPoint_[ie.cellLabel()]);
    }

    label currentDualCell0 = findDualCell(ie.cellLabel(), e[0]);
    label currentDualCell1 = findDualCell(ie.cellLabel(), e[1]);

    ++ie;

    while (true)
    {
        label facei = ie.faceLabel();

        // Mark face as visited.
        doneEFaces[findIndex(eFaces, facei)] = true;

        if (faceToDualPoint_[facei] != -1)
        {
            verts.append(faceToDualPoint_[facei]);
        }

        label celli = ie.cellLabel();

        if (celli == -1)
        {
            // At ending boundary face. We've stored the face point above
            // so this is the whole face.
            break;
        }


        label dualCell0 = findDualCell(celli, e[0]);
        label dualCell1 = findDualCell(celli, e[1]);

        // Generate face. (always if splitFace=true; only if needed to
        // separate cells otherwise)
        if
        (
            splitFace
         || (
                dualCell0 != currentDualCell0
             || dualCell1 != currentDualCell1
            )
        )
        {
            // Close current face.
            addInternalFace
            (
                -1,         // masterPointi
                edgeI,      // masterEdgeI
                -1,         // masterFacei
                edgeOrder,
                currentDualCell0,
                currentDualCell1,
                verts.shrink(),
                meshMod
            );

            // Restart
            currentDualCell0 = dualCell0;
            currentDualCell1 = dualCell1;

            verts.clear();
            if (edgeToDualPoint_[edgeI] != -1)
            {
                verts.append(edgeToDualPoint_[edgeI]);
            }
            if (faceToDualPoint_[facei] != -1)
            {
                verts.append(faceToDualPoint_[facei]);
            }
        }

        if (cellToDualPoint_[celli] != -1)
        {
            verts.append(cellToDualPoint_[celli]);
        }

        ++ie;

        if (ie == ie.end())
        {
            // Back at start face (for internal edge only). See if this needs
            // adding.
            if (isBoundaryEdge.get(edgeI) == 0)
            {
                label startDual = faceToDualPoint_[startFaceLabel];

                if (startDual != -1 && findIndex(verts, startDual) == -1)
                {
                    verts.append(startDual);
                }
            }
            break;
        }
    }

    verts.shrink();
    addInternalFace
    (
        -1,         // masterPointi
        edgeI,      // masterEdgeI
        -1,         // masterFacei
        edgeOrder,
        currentDualCell0,
        currentDualCell1,
        verts,
        meshMod
    );
}


// Walks around circumference of facei. Creates single face. Gets given
// starting (feature) edge to start from. Returns ending edge. (all edges
// in form of index in faceEdges)
void Foam::meshDualiser::createFaceFromInternalFace
(
    const label facei,
    label& fp,
    polyTopoChange& meshMod
) const
{
    const face& f = mesh_.faces()[facei];
    const labelList& fEdges = mesh_.faceEdges()[facei];
    label own = mesh_.faceOwner()[facei];
    label nei = mesh_.faceNeighbour()[facei];

    //Pout<< "createFaceFromInternalFace : At face:" << facei
    //    << " verts:" << f
    //    << " points:" << UIndirectList<point>(mesh_.points(), f)()
    //    << " started walking at edge:" << fEdges[fp]
    //    << " verts:" << mesh_.edges()[fEdges[fp]]
    //    << endl;


    // Walk and collect face.
    DynamicList<label> verts(100);

    verts.append(faceToDualPoint_[facei]);
    verts.append(edgeToDualPoint_[fEdges[fp]]);

    // Step to vertex after edge mid
    fp = f.fcIndex(fp);

    // Get cells on either side of face at that point
    label currentDualCell0 = findDualCell(own, f[fp]);
    label currentDualCell1 = findDualCell(nei, f[fp]);

    forAll(f, i)
    {
        // Check vertex
        if (pointToDualPoint_[f[fp]] != -1)
        {
            verts.append(pointToDualPoint_[f[fp]]);
        }

        // Edge between fp and fp+1
        label edgeI = fEdges[fp];

        if (edgeToDualPoint_[edgeI] != -1)
        {
            verts.append(edgeToDualPoint_[edgeI]);
        }

        // Next vertex on edge
        label nextFp = f.fcIndex(fp);

        // Get dual cells on nextFp to check whether face needs closing.
        label dualCell0 = findDualCell(own, f[nextFp]);
        label dualCell1 = findDualCell(nei, f[nextFp]);

        if (dualCell0 != currentDualCell0 || dualCell1 != currentDualCell1)
        {
            // Check: make sure that there is a midpoint on the edge.
            if (edgeToDualPoint_[edgeI] == -1)
            {
                FatalErrorInFunction
                    << "face:" << facei << " verts:" << f
                    << " points:" << UIndirectList<point>(mesh_.points(), f)()
                    << " no feature edge between " << f[fp]
                    << " and " << f[nextFp] << " although have different"
                    << " dual cells." << endl
                    << "point " << f[fp] << " has dual cells "
                    << currentDualCell0 << " and " << currentDualCell1
                    << " ; point "<< f[nextFp] << " has dual cells "
                    << dualCell0 << " and " << dualCell1
                    << abort(FatalError);
            }


            // Close current face.
            verts.shrink();
            addInternalFace
            (
                -1,         // masterPointi
                -1,         // masterEdgeI
                facei,      // masterFacei
                true,       // edgeOrder,
                currentDualCell0,
                currentDualCell1,
                verts,
                meshMod
            );
            break;
        }

        fp = nextFp;
    }
}


// Given a point on a face converts the faces around the point.
// (pointFaces()). Gets starting face and marks off visited faces in donePFaces.
void Foam::meshDualiser::createFacesAroundBoundaryPoint
(
    const label patchi,
    const label patchPointi,
    const label startFacei,
    polyTopoChange& meshMod,
    boolList& donePFaces            // pFaces visited
) const
{
    const polyBoundaryMesh& patches = mesh_.boundaryMesh();
    const polyPatch& pp = patches[patchi];
    const labelList& pFaces = pp.pointFaces()[patchPointi];
    const labelList& own = mesh_.faceOwner();

    label pointi = pp.meshPoints()[patchPointi];

    if (pointToDualPoint_[pointi] == -1)
    {
        // Not a feature point. Loop over all connected
        // pointFaces.

        // Starting face
        label facei = startFacei;

        DynamicList<label> verts(4);

        while (true)
        {
            label index = findIndex(pFaces, facei-pp.start());

            // Has face been visited already?
            if (donePFaces[index])
            {
                break;
            }
            donePFaces[index] = true;

            // Insert face centre
            verts.append(faceToDualPoint_[facei]);

            label dualCelli = findDualCell(own[facei], pointi);

            // Get the edge before the patchPointi
            const face& f = mesh_.faces()[facei];
            label fp = findIndex(f, pointi);
            label prevFp = f.rcIndex(fp);
            label edgeI = mesh_.faceEdges()[facei][prevFp];

            if (edgeToDualPoint_[edgeI] != -1)
            {
                verts.append(edgeToDualPoint_[edgeI]);
            }

            // Get next boundary face (whilst staying on edge).
            edgeFaceCirculator circ
            (
                mesh_,
                facei,
                true,   // ownerSide
                prevFp, // index of edge in face
                true    // isBoundaryEdge
            );

            do
            {
                ++circ;
            }
            while (mesh_.isInternalFace(circ.faceLabel()));

            // Step to next face
            facei = circ.faceLabel();

            if (facei < pp.start() || facei >= pp.start()+pp.size())
            {
                FatalErrorInFunction
                    << "Walked from face on patch:" << patchi
                    << " to face:" << facei
                    << " fc:" << mesh_.faceCentres()[facei]
                    << " on patch:" << patches.whichPatch(facei)
                    << abort(FatalError);
            }

            // Check if different cell.
            if (dualCelli != findDualCell(own[facei], pointi))
            {
                FatalErrorInFunction
                    << "Different dual cells but no feature edge"
                    << " inbetween point:" << pointi
                    << " coord:" << mesh_.points()[pointi]
                    << abort(FatalError);
            }
        }

        verts.shrink();

        label dualCelli = findDualCell(own[facei], pointi);

        //Bit dodgy: create dualface from the last face (instead of from
        // the central point). This will also use the original faceZone to
        // put the new face (which might span multiple original faces) in.

        addBoundaryFace
        (
            //pointi,     // masterPointi
            -1,         // masterPointi
            -1,         // masterEdgeI
            facei,      // masterFacei
            dualCelli,
            patchi,
            verts,
            meshMod
        );
    }
    else
    {
        label facei = startFacei;

        // Storage for face
        DynamicList<label> verts(mesh_.faces()[facei].size());

        // Starting point.
        verts.append(pointToDualPoint_[pointi]);

        // Find edge between pointi and next point on face.
        const labelList& fEdges = mesh_.faceEdges()[facei];
        label nextEdgeI = fEdges[findIndex(mesh_.faces()[facei], pointi)];
        if (edgeToDualPoint_[nextEdgeI] != -1)
        {
            verts.append(edgeToDualPoint_[nextEdgeI]);
        }

        do
        {
            label index = findIndex(pFaces, facei-pp.start());

            // Has face been visited already?
            if (donePFaces[index])
            {
                break;
            }
            donePFaces[index] = true;

            // Face centre
            verts.append(faceToDualPoint_[facei]);

            // Find edge before pointi on facei
            const labelList& fEdges = mesh_.faceEdges()[facei];
            const face& f = mesh_.faces()[facei];
            label prevFp = f.rcIndex(findIndex(f, pointi));
            label edgeI = fEdges[prevFp];

            if (edgeToDualPoint_[edgeI] != -1)
            {
                // Feature edge. Close any face so far. Note: uses face to
                // create dualFace from. Could use pointi instead.
                verts.append(edgeToDualPoint_[edgeI]);
                addBoundaryFace
                (
                    -1,     // masterPointi
                    -1,     // masterEdgeI
                    facei,  // masterFacei
                    findDualCell(own[facei], pointi),
                    patchi,
                    verts.shrink(),
                    meshMod
                );
                verts.clear();

                verts.append(pointToDualPoint_[pointi]);
                verts.append(edgeToDualPoint_[edgeI]);
            }

            // Cross edgeI to next boundary face
            edgeFaceCirculator circ
            (
                mesh_,
                facei,
                true,   // ownerSide
                prevFp, // index of edge in face
                true    // isBoundaryEdge
            );

            do
            {
                ++circ;
            }
            while (mesh_.isInternalFace(circ.faceLabel()));

            // Step to next face. Quit if not on same patch.
            facei = circ.faceLabel();
        }
        while
        (
            facei != startFacei
         && facei >= pp.start()
         && facei < pp.start()+pp.size()
        );

        if (verts.size() > 2)
        {
            // Note: face created from face, not from pointi
            addBoundaryFace
            (
                -1,             // masterPointi
                -1,             // masterEdgeI
                startFacei,     // masterFacei
                findDualCell(own[facei], pointi),
                patchi,
                verts.shrink(),
                meshMod
            );
        }
    }
}


// * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //

Foam::meshDualiser::meshDualiser(const polyMesh& mesh)
:
    mesh_(mesh),
    pointToDualCells_(mesh_.nPoints()),
    pointToDualPoint_(mesh_.nPoints(), -1),
    cellToDualPoint_(mesh_.nCells()),
    faceToDualPoint_(mesh_.nFaces(), -1),
    edgeToDualPoint_(mesh_.nEdges(), -1)
{}


// * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //

void Foam::meshDualiser::setRefinement
(
    const bool splitFace,
    const labelList& featureFaces,
    const labelList& featureEdges,
    const labelList& singleCellFeaturePoints,
    const labelList& multiCellFeaturePoints,
    polyTopoChange& meshMod
)
{
    const labelList& own = mesh_.faceOwner();
    const labelList& nei = mesh_.faceNeighbour();
    const vectorField& cellCentres = mesh_.cellCentres();

    // Mark boundary edges and points.
    // (Note: in 1.4.2 we can use the built-in mesh point ordering
    //  facility instead)
    PackedBoolList isBoundaryEdge(mesh_.nEdges());
    for (label facei = mesh_.nInternalFaces(); facei < mesh_.nFaces(); facei++)
    {
        const labelList& fEdges = mesh_.faceEdges()[facei];

        forAll(fEdges, i)
        {
            isBoundaryEdge.set(fEdges[i], 1);
        }
    }


    if (splitFace)
    {
        // This is a special mode where whenever we are walking around an edge
        // every area through a cell becomes a separate dualface. So two
        // dual cells will probably have more than one dualface between them!
        // This mode implies that
        // - all faces have to be feature faces since there has to be a
        //   dualpoint at the face centre.
        // - all edges have to be feature edges ,,
        boolList featureFaceSet(mesh_.nFaces(), false);
        forAll(featureFaces, i)
        {
            featureFaceSet[featureFaces[i]] = true;
        }
        label facei = findIndex(featureFaceSet, false);

        if (facei != -1)
        {
            FatalErrorInFunction
                << "In split-face-mode (splitFace=true) but not all faces"
                << " marked as feature faces." << endl
                << "First conflicting face:" << facei
                << " centre:" << mesh_.faceCentres()[facei]
                << abort(FatalError);
        }

        boolList featureEdgeSet(mesh_.nEdges(), false);
        forAll(featureEdges, i)
        {
            featureEdgeSet[featureEdges[i]] = true;
        }
        label edgeI = findIndex(featureEdgeSet, false);

        if (edgeI != -1)
        {
            const edge& e = mesh_.edges()[edgeI];
            FatalErrorInFunction
                << "In split-face-mode (splitFace=true) but not all edges"
                << " marked as feature edges." << endl
                << "First conflicting edge:" << edgeI
                << " verts:" << e
                << " coords:" << mesh_.points()[e[0]] << mesh_.points()[e[1]]
                << abort(FatalError);
        }
    }
    else
    {
        // Check that all boundary faces are feature faces.

        boolList featureFaceSet(mesh_.nFaces(), false);
        forAll(featureFaces, i)
        {
            featureFaceSet[featureFaces[i]] = true;
        }
        for
        (
            label facei = mesh_.nInternalFaces();
            facei < mesh_.nFaces();
            facei++
        )
        {
            if (!featureFaceSet[facei])
            {
                FatalErrorInFunction
                    << "Not all boundary faces marked as feature faces."
                    << endl
                    << "First conflicting face:" << facei
                    << " centre:" << mesh_.faceCentres()[facei]
                    << abort(FatalError);
            }
        }
    }




    // Start creating cells, points, and faces (in that order)


    // 1. Mark which cells to create
    // Mostly every point becomes one cell but sometimes (for feature points)
    // all cells surrounding a feature point become cells. Also a non-manifold
    // point can create two cells! So a dual cell is uniquely defined by a
    // mesh point + cell (as in pointCells index)

    // 2. Mark which face centres to create

    // 3. Internal faces can now consist of
    //      - only cell centres of walk around edge
    //      - cell centres + face centres of walk around edge
    //      - same but now other side is not a single cell

    // 4. Boundary faces (or internal faces between cell zones!) now consist of
    //      - walk around boundary point.



    autoPtr<OFstream> dualCcStr;
    if (debug)
    {
        dualCcStr.reset(new OFstream("dualCc.obj"));
        Pout<< "Dumping centres of dual cells to " << dualCcStr().name()
            << endl;
    }


    // Dual cells (from points)
    // ~~~~~~~~~~~~~~~~~~~~~~~~

    // pointToDualCells_[pointi]
    // - single entry : all cells surrounding point all become the same
    //                  cell.
    // - multiple entries: in order of pointCells.


    // feature points that become single cell
    forAll(singleCellFeaturePoints, i)
    {
        label pointi = singleCellFeaturePoints[i];

        pointToDualPoint_[pointi] = meshMod.addPoint
        (
            mesh_.points()[pointi],
            pointi,                                 // masterPoint
            mesh_.pointZones().whichZone(pointi),   // zoneID
            true                                    // inCell
        );

        // Generate single cell
        pointToDualCells_[pointi].setSize(1);
        pointToDualCells_[pointi][0] = meshMod.addCell
        (
            pointi, //masterPointID,
            -1,     //masterEdgeID,
            -1,     //masterFaceID,
            -1,     //masterCellID,
            -1      //zoneID
        );
        if (dualCcStr.valid())
        {
            meshTools::writeOBJ(dualCcStr(), mesh_.points()[pointi]);
        }
    }

    // feature points that become multiple cells
    forAll(multiCellFeaturePoints, i)
    {
        label pointi = multiCellFeaturePoints[i];

        if (pointToDualCells_[pointi].size() > 0)
        {
            FatalErrorInFunction
                << "Point " << pointi << " at:" << mesh_.points()[pointi]
                << " is both in singleCellFeaturePoints"
                << " and multiCellFeaturePoints."
                << abort(FatalError);
        }

        pointToDualPoint_[pointi] = meshMod.addPoint
        (
            mesh_.points()[pointi],
            pointi,                                 // masterPoint
            mesh_.pointZones().whichZone(pointi),   // zoneID
            true                                    // inCell
        );

        // Create dualcell for every cell connected to dual point

        const labelList& pCells = mesh_.pointCells()[pointi];

        pointToDualCells_[pointi].setSize(pCells.size());

        forAll(pCells, pCelli)
        {
            pointToDualCells_[pointi][pCelli] = meshMod.addCell
            (
                pointi,                                     //masterPointID
                -1,                                         //masterEdgeID
                -1,                                         //masterFaceID
                -1,                                         //masterCellID
                mesh_.cellZones().whichZone(pCells[pCelli]) //zoneID
            );
            if (dualCcStr.valid())
            {
                meshTools::writeOBJ
                (
                    dualCcStr(),
                    0.5*(mesh_.points()[pointi]+cellCentres[pCells[pCelli]])
                );
            }
        }
    }
    // Normal points
    forAll(mesh_.points(), pointi)
    {
        if (pointToDualCells_[pointi].empty())
        {
            pointToDualCells_[pointi].setSize(1);
            pointToDualCells_[pointi][0] = meshMod.addCell
            (
                pointi, //masterPointID,
                -1,     //masterEdgeID,
                -1,     //masterFaceID,
                -1,     //masterCellID,
                -1      //zoneID
            );

            if (dualCcStr.valid())
            {
                meshTools::writeOBJ(dualCcStr(), mesh_.points()[pointi]);
            }
        }
    }


    // Dual points (from cell centres, feature faces, feature edges)
    // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~

    forAll(cellToDualPoint_, celli)
    {
        cellToDualPoint_[celli] = meshMod.addPoint
        (
            cellCentres[celli],
            mesh_.faces()[mesh_.cells()[celli][0]][0],     // masterPoint
            -1,     // zoneID
            true    // inCell
        );
    }

    // From face to dual point

    forAll(featureFaces, i)
    {
        label facei = featureFaces[i];

        faceToDualPoint_[facei] = meshMod.addPoint
        (
            mesh_.faceCentres()[facei],
            mesh_.faces()[facei][0],     // masterPoint
            -1,     // zoneID
            true    // inCell
        );
    }
    // Detect whether different dual cells on either side of a face. This
    // would neccesitate having a dual face built from the face and thus a
    // dual point at the face centre.
    for (label facei = 0; facei < mesh_.nInternalFaces(); facei++)
    {
        if (faceToDualPoint_[facei] == -1)
        {
            const face& f = mesh_.faces()[facei];

            forAll(f, fp)
            {
                label ownDualCell = findDualCell(own[facei], f[fp]);
                label neiDualCell = findDualCell(nei[facei], f[fp]);

                if (ownDualCell != neiDualCell)
                {
                    faceToDualPoint_[facei] = meshMod.addPoint
                    (
                        mesh_.faceCentres()[facei],
                        f[fp],  // masterPoint
                        -1,     // zoneID
                        true    // inCell
                    );

                    break;
                }
            }
        }
    }

    // From edge to dual point

    forAll(featureEdges, i)
    {
        label edgeI = featureEdges[i];

        const edge& e = mesh_.edges()[edgeI];

        edgeToDualPoint_[edgeI] = meshMod.addPoint
        (
            e.centre(mesh_.points()),
            e[0],   // masterPoint
            -1,     // zoneID
            true    // inCell
        );
    }

    // Detect whether different dual cells on either side of an edge. This
    // would neccesitate having a dual face built perpendicular to the edge
    // and thus a dual point at the mid of the edge.
    // Note: not really true - the face can be built without the edge centre!
    const labelListList& edgeCells = mesh_.edgeCells();

    forAll(edgeCells, edgeI)
    {
       if (edgeToDualPoint_[edgeI] == -1)
       {
            const edge& e = mesh_.edges()[edgeI];

            // We need a point on the edge if not all cells on both sides
            // are the same.

            const labelList& eCells = mesh_.edgeCells()[edgeI];

            label dualE0 = findDualCell(eCells[0], e[0]);
            label dualE1 = findDualCell(eCells[0], e[1]);

            for (label i = 1; i < eCells.size(); i++)
            {
                label newDualE0 = findDualCell(eCells[i], e[0]);

                if (dualE0 != newDualE0)
                {
                    edgeToDualPoint_[edgeI] = meshMod.addPoint
                    (
                        e.centre(mesh_.points()),
                        e[0],                               // masterPoint
                        mesh_.pointZones().whichZone(e[0]), // zoneID
                        true                                // inCell
                    );

                    break;
                }

                label newDualE1 = findDualCell(eCells[i], e[1]);

                if (dualE1 != newDualE1)
                {
                    edgeToDualPoint_[edgeI] = meshMod.addPoint
                    (
                        e.centre(mesh_.points()),
                        e[1],   // masterPoint
                        mesh_.pointZones().whichZone(e[1]), // zoneID
                        true    // inCell
                    );

                    break;
                }
            }
        }
    }

    // Make sure all boundary edges emanating from feature points are
    // feature edges as well.
    forAll(singleCellFeaturePoints, i)
    {
        generateDualBoundaryEdges
        (
            isBoundaryEdge,
            singleCellFeaturePoints[i],
            meshMod
        );
    }
    forAll(multiCellFeaturePoints, i)
    {
        generateDualBoundaryEdges
        (
            isBoundaryEdge,
            multiCellFeaturePoints[i],
            meshMod
        );
    }


    // Check for duplicate points
    if (debug)
    {
        dumpPolyTopoChange(meshMod, "generatedPoints_");
        checkPolyTopoChange(meshMod);
    }


    // Now we have all points and cells
    // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
    //  - pointToDualCells_ : per point a single dualCell or multiple dualCells
    //  - pointToDualPoint_ : per point -1 or the dual point at the coordinate
    //  - edgeToDualPoint_  : per edge -1 or the edge centre
    //  - faceToDualPoint_  : per face -1 or the face centre
    //  - cellToDualPoint_  : per cell the cell centre
    // Now we have to walk all edges and construct faces. Either single face
    // per edge or multiple (-if nonmanifold edge -if different dualcells)

    const edgeList& edges = mesh_.edges();

    forAll(edges, edgeI)
    {
        const labelList& eFaces = mesh_.edgeFaces()[edgeI];

        boolList doneEFaces(eFaces.size(), false);

        forAll(eFaces, i)
        {
            if (!doneEFaces[i])
            {
                // We found a face that hasn't yet been visited. This might
                // happen for non-manifold edges where a single edge can
                // become multiple faces.

                label startFacei = eFaces[i];

                //Pout<< "Walking edge:" << edgeI
                //    << " points:" << mesh_.points()[e[0]]
                //    << mesh_.points()[e[1]]
                //    << " startFace:" << startFacei
                //    << " at:" << mesh_.faceCentres()[startFacei]
                //    << endl;

                createFacesAroundEdge
                (
                    splitFace,
                    isBoundaryEdge,
                    edgeI,
                    startFacei,
                    meshMod,
                    doneEFaces
                );
            }
        }
    }

    if (debug)
    {
        dumpPolyTopoChange(meshMod, "generatedFacesFromEdges_");
    }

    // Create faces from feature faces. These can be internal or external faces.
    // - feature face : centre needs to be included.
    //      - single cells on either side: triangulate
    //      - multiple cells: create single face between unique cell pair. Only
    //                        create face where cells differ on either side.
    // - non-feature face : inbetween cell zones.
    forAll(faceToDualPoint_, facei)
    {
        if (faceToDualPoint_[facei] != -1 && mesh_.isInternalFace(facei))
        {
            const face& f = mesh_.faces()[facei];
            const labelList& fEdges = mesh_.faceEdges()[facei];

            // Starting edge
            label fp = 0;

            do
            {
                // Find edge that is in dual mesh and where the
                // next point (fp+1) has different dual cells on either side.
                bool foundStart = false;

                do
                {
                    if
                    (
                        edgeToDualPoint_[fEdges[fp]] != -1
                    && !sameDualCell(facei, f.nextLabel(fp))
                    )
                    {
                        foundStart = true;
                        break;
                    }
                    fp = f.fcIndex(fp);
                }
                while (fp != 0);

                if (!foundStart)
                {
                    break;
                }

                // Walk from edge fp and generate a face.
                createFaceFromInternalFace
                (
                    facei,
                    fp,
                    meshMod
                );
            }
            while (fp != 0);
        }
    }

    if (debug)
    {
        dumpPolyTopoChange(meshMod, "generatedFacesFromFeatFaces_");
    }


    // Create boundary faces. Every boundary point has one or more dualcells.
    // These need to be closed.
    const polyBoundaryMesh& patches = mesh_.boundaryMesh();

    forAll(patches, patchi)
    {
        const polyPatch& pp = patches[patchi];

        const labelListList& pointFaces = pp.pointFaces();

        forAll(pointFaces, patchPointi)
        {
            const labelList& pFaces = pointFaces[patchPointi];

            boolList donePFaces(pFaces.size(), false);

            forAll(pFaces, i)
            {
                if (!donePFaces[i])
                {
                    // Starting face
                    label startFacei = pp.start()+pFaces[i];

                    //Pout<< "Walking around point:" << pointi
                    //    << " coord:" << mesh_.points()[pointi]
                    //    << " on patch:" << patchi
                    //    << " startFace:" << startFacei
                    //    << " at:" << mesh_.faceCentres()[startFacei]
                    //    << endl;

                    createFacesAroundBoundaryPoint
                    (
                        patchi,
                        patchPointi,
                        startFacei,
                        meshMod,
                        donePFaces            // pFaces visited
                    );
                }
            }
        }
    }

    if (debug)
    {
        dumpPolyTopoChange(meshMod, "generatedFacesFromBndFaces_");
    }
}


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