File: mergeOrSplitBaffles.C

package info (click to toggle)
openfoam 1812%2Bdfsg1-2
  • links: PTS, VCS
  • area: main
  • in suites: buster
  • size: 220,284 kB
  • sloc: cpp: 1,038,902; sh: 14,536; ansic: 8,240; lex: 657; xml: 387; python: 300; awk: 212; makefile: 94; sed: 88; csh: 3
file content (583 lines) | stat: -rw-r--r-- 16,626 bytes parent folder | download
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
   \\    /   O peration     |
    \\  /    A nd           | Copyright (C) 2011-2016 OpenFOAM Foundation
     \\/     M anipulation  | Copyright (C) 2016-2018 OpenCFD Ltd.
-------------------------------------------------------------------------------
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/>.

Application
    mergeOrSplitBaffles

Group
    grpMeshManipulationUtilities

Description
    Detects boundary faces that share points (baffles). Either merges them or
    duplicate the points.

Usage
    \b mergeOrSplitBaffles [OPTION]

    Options:
      - \par -detectOnly
        Detect baffles and write to faceSet duplicateFaces.

      - \par -split
        Detect baffles and duplicate the points (used so the two sides
        can move independently)

      - \par -dict \<dictionary\>
        Specify a dictionary to read actions from.

Note
    - can only handle pairwise boundary faces. So three faces using
      the same points is not handled (is illegal mesh anyway)

    - surfaces consisting of duplicate faces can be topologically split
    if the points on the interior of the surface cannot walk to all the
    cells that use them in one go.

    - Parallel operation (where duplicate face is perpendicular to a coupled
    boundary) is supported but not really tested.
    (Note that coupled faces themselves are not seen as duplicate faces)

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

#include "argList.H"
#include "Time.H"
#include "syncTools.H"
#include "faceSet.H"
#include "pointSet.H"
#include "meshTools.H"
#include "polyTopoChange.H"
#include "polyRemoveFace.H"
#include "polyModifyFace.H"
#include "indirectPrimitivePatch.H"
#include "processorPolyPatch.H"
#include "localPointRegion.H"
#include "duplicatePoints.H"
#include "ReadFields.H"
#include "volFields.H"
#include "surfaceFields.H"
#include "processorMeshes.H"

using namespace Foam;

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

void insertDuplicateMerge
(
    const polyMesh& mesh,
    const labelList& boundaryFaces,
    const labelList& duplicates,
    polyTopoChange& meshMod
)
{
    const faceList& faces = mesh.faces();
    const labelList& faceOwner = mesh.faceOwner();
    const faceZoneMesh& faceZones = mesh.faceZones();

    forAll(duplicates, bFacei)
    {
        label otherFacei = duplicates[bFacei];

        if (otherFacei != -1 && otherFacei > bFacei)
        {
            // Two duplicate faces. Merge.

            label face0 = boundaryFaces[bFacei];
            label face1 = boundaryFaces[otherFacei];

            label own0 = faceOwner[face0];
            label own1 = faceOwner[face1];

            if (own0 < own1)
            {
                // Use face0 as the new internal face.
                label zoneID = faceZones.whichZone(face0);
                bool zoneFlip = false;

                if (zoneID >= 0)
                {
                    const faceZone& fZone = faceZones[zoneID];
                    zoneFlip = fZone.flipMap()[fZone.whichFace(face0)];
                }

                meshMod.setAction(polyRemoveFace(face1));
                meshMod.setAction
                (
                    polyModifyFace
                    (
                        faces[face0],           // modified face
                        face0,                  // label of face being modified
                        own0,                   // owner
                        own1,                   // neighbour
                        false,                  // face flip
                        -1,                     // patch for face
                        false,                  // remove from zone
                        zoneID,                 // zone for face
                        zoneFlip                // face flip in zone
                    )
                );
            }
            else
            {
                // Use face1 as the new internal face.
                label zoneID = faceZones.whichZone(face1);
                bool zoneFlip = false;

                if (zoneID >= 0)
                {
                    const faceZone& fZone = faceZones[zoneID];
                    zoneFlip = fZone.flipMap()[fZone.whichFace(face1)];
                }

                meshMod.setAction(polyRemoveFace(face0));
                meshMod.setAction
                (
                    polyModifyFace
                    (
                        faces[face1],           // modified face
                        face1,                  // label of face being modified
                        own1,                   // owner
                        own0,                   // neighbour
                        false,                  // face flip
                        -1,                     // patch for face
                        false,                  // remove from zone
                        zoneID,                 // zone for face
                        zoneFlip                // face flip in zone
                    )
                );
            }
        }
    }
}


label patchSize(const polyMesh& mesh, const labelList& patchIDs)
{
    const polyBoundaryMesh& patches = mesh.boundaryMesh();

    label sz = 0;
    forAll(patchIDs, i)
    {
        const polyPatch& pp = patches[patchIDs[i]];
        sz += pp.size();
    }
    return sz;
}


labelList patchFaces(const polyMesh& mesh, const labelList& patchIDs)
{
    const polyBoundaryMesh& patches = mesh.boundaryMesh();

    labelList faceIDs(patchSize(mesh, patchIDs));
    label sz = 0;
    forAll(patchIDs, i)
    {
        const polyPatch& pp = patches[patchIDs[i]];

        forAll(pp, ppi)
        {
            faceIDs[sz++] = pp.start()+ppi;
        }
    }

    if (faceIDs.size() != sz)
    {
        FatalErrorInFunction << exit(FatalError);
    }

    return faceIDs;
}


labelList findBaffles(const polyMesh& mesh, const labelList& boundaryFaces)
{
    // Get all duplicate face labels (in boundaryFaces indices!).
    labelList duplicates = localPointRegion::findDuplicateFaces
    (
        mesh,
        boundaryFaces
    );


    // Check that none are on processor patches
    const polyBoundaryMesh& patches = mesh.boundaryMesh();

    forAll(duplicates, bFacei)
    {
        if (duplicates[bFacei] != -1)
        {
            label facei = boundaryFaces[bFacei];
            label patchi = patches.whichPatch(facei);

            if (isA<processorPolyPatch>(patches[patchi]))
            {
                FatalErrorInFunction
                    << "Duplicate face " << facei
                    << " is on a processorPolyPatch."
                    << "This is not allowed." << nl
                    << "Face:" << facei
                    << " is on patch:" << patches[patchi].name()
                    << abort(FatalError);
            }
        }
    }


    // Write to faceSet for ease of post-processing.
    {
        faceSet duplicateSet
        (
            mesh,
            "duplicateFaces",
            mesh.nBoundaryFaces()/256
        );

        forAll(duplicates, bFacei)
        {
            label otherFacei = duplicates[bFacei];

            if (otherFacei != -1 && otherFacei > bFacei)
            {
                duplicateSet.insert(boundaryFaces[bFacei]);
                duplicateSet.insert(boundaryFaces[otherFacei]);
            }
        }

        Info<< "Writing " << returnReduce(duplicateSet.size(), sumOp<label>())
            << " duplicate faces to faceSet " << duplicateSet.objectPath()
            << nl << endl;
        duplicateSet.write();
    }

    return duplicates;
}


int main(int argc, char *argv[])
{
    argList::addNote
    (
        "Detect faces that share points (baffles).\n"
        "Merge them or duplicate the points."
    );

    #include "addOverwriteOption.H"
    #include "addRegionOption.H"
    #include "addDictOption.H"
    argList::addBoolOption
    (
        "detectOnly",
        "Find baffles only, but do not merge or split them"
    );
    argList::addBoolOption
    (
        "split",
        "Topologically split duplicate surfaces"
    );

    argList::noFunctionObjects();  // Never use function objects

    #include "setRootCase.H"
    #include "createTime.H"
    #include "createNamedMesh.H"

    const word oldInstance = mesh.pointsInstance();
    const polyBoundaryMesh& patches = mesh.boundaryMesh();

    const bool readDict   = args.found("dict");
    const bool split      = args.found("split");
    const bool overwrite  = args.found("overwrite");
    const bool detectOnly = args.found("detectOnly");

    if (readDict && (split || detectOnly))
    {
        FatalErrorInFunction
            << "Use of dictionary for settings not compatible with"
            << " using command line arguments for \"split\""
            << " or \"detectOnly\"" << exit(FatalError);
    }


    labelList detectPatchIDs;
    labelList splitPatchIDs;
    labelList mergePatchIDs;

    if (readDict)
    {
        const word dictName;
        #include "setSystemMeshDictionaryIO.H"

        Info<< "Reading " << dictName << "\n" << endl;
        IOdictionary dict(dictIO);

        if (dict.found("detect"))
        {
            detectPatchIDs = patches.patchSet
            (
                dict.subDict("detect").get<wordRes>("patches")
            ).sortedToc();

            Info<< "Detecting baffles on " << detectPatchIDs.size()
                << " patches with "
                << returnReduce(patchSize(mesh, detectPatchIDs), sumOp<label>())
                << " faces" << endl;
        }
        if (dict.found("merge"))
        {
            mergePatchIDs = patches.patchSet
            (
                dict.subDict("merge").get<wordRes>("patches")
            ).sortedToc();

            Info<< "Detecting baffles on " << mergePatchIDs.size()
                << " patches with "
                << returnReduce(patchSize(mesh, mergePatchIDs), sumOp<label>())
                << " faces" << endl;
        }
        if (dict.found("split"))
        {
            splitPatchIDs = patches.patchSet
            (
                dict.subDict("split").get<wordRes>("patches")
            ).sortedToc();

            Info<< "Detecting baffles on " << splitPatchIDs.size()
                << " patches with "
                << returnReduce(patchSize(mesh, splitPatchIDs), sumOp<label>())
                << " faces" << endl;
        }
    }
    else
    {
        if (detectOnly)
        {
            detectPatchIDs = identity(patches.size());
        }
        else if (split)
        {
            splitPatchIDs = identity(patches.size());
        }
        else
        {
            mergePatchIDs = identity(patches.size());
        }
    }


    if (detectPatchIDs.size())
    {
        findBaffles(mesh, patchFaces(mesh, detectPatchIDs));

        if (detectOnly)
        {
            return 0;
        }
    }



    // Read objects in time directory
    IOobjectList objects(mesh, runTime.timeName());

    // Read vol fields.

    PtrList<volScalarField> vsFlds;
    ReadFields(mesh, objects, vsFlds);

    PtrList<volVectorField> vvFlds;
    ReadFields(mesh, objects, vvFlds);

    PtrList<volSphericalTensorField> vstFlds;
    ReadFields(mesh, objects, vstFlds);

    PtrList<volSymmTensorField> vsymtFlds;
    ReadFields(mesh, objects, vsymtFlds);

    PtrList<volTensorField> vtFlds;
    ReadFields(mesh, objects, vtFlds);

    // Read surface fields.

    PtrList<surfaceScalarField> ssFlds;
    ReadFields(mesh, objects, ssFlds);

    PtrList<surfaceVectorField> svFlds;
    ReadFields(mesh, objects, svFlds);

    PtrList<surfaceSphericalTensorField> sstFlds;
    ReadFields(mesh, objects, sstFlds);

    PtrList<surfaceSymmTensorField> ssymtFlds;
    ReadFields(mesh, objects, ssymtFlds);

    PtrList<surfaceTensorField> stFlds;
    ReadFields(mesh, objects, stFlds);



    if (mergePatchIDs.size())
    {
        Info<< "Merging duplicate faces" << nl << endl;

        // Mesh change engine
        polyTopoChange meshMod(mesh);

        const labelList boundaryFaces(patchFaces(mesh, mergePatchIDs));

        // Get all duplicate face pairs (in boundaryFaces indices!).
        labelList duplicates(findBaffles(mesh, boundaryFaces));

        // Merge into internal faces.
        insertDuplicateMerge(mesh, boundaryFaces, duplicates, meshMod);

        if (!overwrite)
        {
            ++runTime;
        }

        // Change the mesh. No inflation.
        autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);

        // Update fields
        mesh.updateMesh(map());

        // Move mesh (since morphing does not do this)
        if (map().hasMotionPoints())
        {
            mesh.movePoints(map().preMotionPoints());
        }

        if (overwrite)
        {
            mesh.setInstance(oldInstance);
        }
        Info<< "Writing mesh to time " << runTime.timeName() << endl;
        mesh.write();
    }


    if (splitPatchIDs.size())
    {
        Info<< "Topologically splitting duplicate surfaces"
            << ", i.e. duplicating points internal to duplicate surfaces"
            << nl << endl;

        // Determine points on split patches
        DynamicList<label> candidates;
        {
            label sz = 0;
            forAll(splitPatchIDs, i)
            {
                sz += patches[splitPatchIDs[i]].nPoints();
            }
            candidates.setCapacity(sz);

            bitSet isCandidate(mesh.nPoints());
            forAll(splitPatchIDs, i)
            {
                const labelList& mp = patches[splitPatchIDs[i]].meshPoints();
                forAll(mp, mpi)
                {
                    label pointi = mp[mpi];
                    if (isCandidate.set(pointi))
                    {
                        candidates.append(pointi);
                    }
                }
            }
        }


        // Analyse which points need to be duplicated
        localPointRegion regionSide(mesh, candidates);

        // Point duplication engine
        duplicatePoints pointDuplicator(mesh);

        // Mesh change engine
        polyTopoChange meshMod(mesh);

        // Insert topo changes
        pointDuplicator.setRefinement(regionSide, meshMod);

        if (!overwrite)
        {
            ++runTime;
        }

        // Change the mesh. No inflation.
        autoPtr<mapPolyMesh> map = meshMod.changeMesh(mesh, false);

        // Update fields
        mesh.updateMesh(map());

        // Move mesh (since morphing does not do this)
        if (map().hasMotionPoints())
        {
            mesh.movePoints(map().preMotionPoints());
        }

        if (overwrite)
        {
            mesh.setInstance(oldInstance);
        }
        Info<< "Writing mesh to time " << runTime.timeName() << endl;
        mesh.write();

        topoSet::removeFiles(mesh);
        processorMeshes::removeFiles(mesh);

        // Dump duplicated points (if any)
        const labelList& pointMap = map().pointMap();

        labelList nDupPerPoint(map().nOldPoints(), 0);

        pointSet dupPoints(mesh, "duplicatedPoints", 100);

        forAll(pointMap, pointi)
        {
            label oldPointi = pointMap[pointi];

            nDupPerPoint[oldPointi]++;

            if (nDupPerPoint[oldPointi] > 1)
            {
                dupPoints.insert(map().reversePointMap()[oldPointi]);
                dupPoints.insert(pointi);
            }
        }

        Info<< "Writing " << returnReduce(dupPoints.size(), sumOp<label>())
            << " duplicated points to pointSet "
            << dupPoints.objectPath() << nl << endl;

        dupPoints.write();
    }

    Info<< "End\n" << endl;

    return 0;
}


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