File: writeMeshObj.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 (580 lines) | stat: -rw-r--r-- 14,086 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
/*---------------------------------------------------------------------------*\
  =========                 |
  \\      /  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/>.

Application
    writeMeshObj

Description
    For mesh debugging: writes mesh as three separate OBJ files which can
    be viewed with e.g. javaview.

    meshPoints_XXX.obj : all points and edges as lines.
    meshFaceCentres_XXX.obj : all face centres.
    meshCellCentres_XXX.obj : all cell centres.

    patch_YYY_XXX.obj : all face centres of patch YYY

    Optional: - patch faces (as polygons) : patchFaces_YYY_XXX.obj
              - non-manifold edges : patchEdges_YYY_XXX.obj

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

#include "argList.H"
#include "timeSelector.H"
#include "Time.H"
#include "polyMesh.H"
#include "OFstream.H"
#include "meshTools.H"
#include "cellSet.H"
#include "faceSet.H"
#include "SubField.H"

using namespace Foam;

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

void writeOBJ(const point& pt, Ostream& os)
{
    os  << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
}

// All edges of mesh
void writePoints(const polyMesh& mesh, const fileName& timeName)
{
    label vertI = 0;

    fileName pointFile(mesh.time().path()/"meshPoints_" + timeName + ".obj");

    Info<< "Writing mesh points and edges to " << pointFile << endl;

    OFstream pointStream(pointFile);

    forAll(mesh.points(), pointi)
    {
        writeOBJ(mesh.points()[pointi], pointStream);
        vertI++;
    }

    forAll(mesh.edges(), edgeI)
    {
        const edge& e = mesh.edges()[edgeI];

        pointStream << "l " << e.start() + 1 << ' ' << e.end() + 1 << nl;
    }
}


// Edges for subset of cells
void writePoints
(
    const polyMesh& mesh,
    const labelList& cellLabels,
    const fileName& timeName
)
{
    fileName fName(mesh.time().path()/"meshPoints_" + timeName + ".obj");

    Info<< "Writing mesh points and edges to " << fName << endl;

    OFstream str(fName);

    // OBJ file vertex
    label vertI = 0;

    // From point to OBJ file vertex
    Map<label> pointToObj(6*cellLabels.size());

    forAll(cellLabels, i)
    {
        label celli = cellLabels[i];

        const labelList& cEdges = mesh.cellEdges()[celli];

        forAll(cEdges, cEdgeI)
        {
            const edge& e = mesh.edges()[cEdges[cEdgeI]];

            label v0;

            Map<label>::iterator e0Fnd = pointToObj.find(e[0]);

            if (e0Fnd == pointToObj.end())
            {
                meshTools::writeOBJ(str, mesh.points()[e[0]]);
                v0 = vertI++;
                pointToObj.insert(e[0], v0);
            }
            else
            {
                v0 = e0Fnd();
            }

            label v1;

            Map<label>::iterator e1Fnd = pointToObj.find(e[1]);

            if (e1Fnd == pointToObj.end())
            {
                meshTools::writeOBJ(str, mesh.points()[e[1]]);
                v1 = vertI++;
                pointToObj.insert(e[1], v1);
            }
            else
            {
                v1 = e1Fnd();
            }


            str << "l " << v0+1 << ' ' << v1+1 << nl;
        }
    }
}


// Edges of single cell
void writePoints
(
    const polyMesh& mesh,
    const label celli,
    const fileName& timeName
)
{
    fileName fName
    (
        mesh.time().path()
      / "meshPoints_" + timeName + '_' + name(celli) + ".obj"
    );

    Info<< "Writing mesh points and edges to " << fName << endl;

    OFstream pointStream(fName);

    const cell& cFaces = mesh.cells()[celli];

    meshTools::writeOBJ(pointStream, mesh.faces(), mesh.points(), cFaces);
}



// All face centres
void writeFaceCentres(const polyMesh& mesh,const fileName& timeName)
{
    fileName faceFile
    (
        mesh.time().path()
      / "meshFaceCentres_" + timeName + ".obj"
    );

    Info<< "Writing mesh face centres to " << faceFile << endl;

    OFstream faceStream(faceFile);

    forAll(mesh.faceCentres(), facei)
    {
        writeOBJ(mesh.faceCentres()[facei], faceStream);
    }
}


void writeCellCentres(const polyMesh& mesh, const fileName& timeName)
{
    fileName cellFile
    (
        mesh.time().path()/"meshCellCentres_" + timeName + ".obj"
    );

    Info<< "Writing mesh cell centres to " << cellFile << endl;

    OFstream cellStream(cellFile);

    forAll(mesh.cellCentres(), celli)
    {
        writeOBJ(mesh.cellCentres()[celli], cellStream);
    }
}


void writePatchCentres
(
    const polyMesh& mesh,
    const fileName& timeName
)
{
    const polyBoundaryMesh& patches = mesh.boundaryMesh();

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

        fileName faceFile
        (
            mesh.time().path()/"patch_" + pp.name() + '_' + timeName + ".obj"
        );

        Info<< "Writing patch face centres to " << faceFile << endl;

        OFstream patchFaceStream(faceFile);

        forAll(pp.faceCentres(), facei)
        {
            writeOBJ(pp.faceCentres()[facei], patchFaceStream);
        }
    }
}


void writePatchFaces
(
    const polyMesh& mesh,
    const fileName& timeName
)
{
    const polyBoundaryMesh& patches = mesh.boundaryMesh();

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

        fileName faceFile
        (
            mesh.time().path()
          / "patchFaces_" + pp.name() + '_' + timeName + ".obj"
        );

        Info<< "Writing patch faces to " << faceFile << endl;

        OFstream patchFaceStream(faceFile);

        forAll(pp.localPoints(), pointi)
        {
            writeOBJ(pp.localPoints()[pointi], patchFaceStream);
        }

        forAll(pp.localFaces(), facei)
        {
            const face& f = pp.localFaces()[facei];

            patchFaceStream<< 'f';

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


void writePatchBoundaryEdges
(
    const polyMesh& mesh,
    const fileName& timeName
)
{
    const polyBoundaryMesh& patches = mesh.boundaryMesh();

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

        fileName edgeFile
        (
            mesh.time().path()
          / "patchEdges_" + pp.name() + '_' + timeName + ".obj"
        );

        Info<< "Writing patch edges to " << edgeFile << endl;

        OFstream patchEdgeStream(edgeFile);

        forAll(pp.localPoints(), pointi)
        {
            writeOBJ(pp.localPoints()[pointi], patchEdgeStream);
        }

        for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
        {
            if (pp.edgeFaces()[edgeI].size() == 1)
            {
                const edge& e = pp.edges()[edgeI];

                patchEdgeStream<< "l " << e[0]+1 << ' ' << e[1]+1 << nl;
            }
        }
    }
}


void writePointCells
(
    const polyMesh& mesh,
    const label pointi,
    const fileName& timeName
)
{
    const labelList& pCells = mesh.pointCells()[pointi];

    labelHashSet allEdges(6*pCells.size());

    forAll(pCells, i)
    {
        const labelList& cEdges = mesh.cellEdges()[pCells[i]];

        forAll(cEdges, i)
        {
            allEdges.insert(cEdges[i]);
        }
    }


    fileName pFile
    (
        mesh.time().path()
      / "pointEdges_" + timeName + '_' + name(pointi) + ".obj"
    );

    Info<< "Writing pointEdges to " << pFile << endl;

    OFstream pointStream(pFile);

    label vertI = 0;

    forAllConstIter(labelHashSet, allEdges, iter)
    {
        const edge& e = mesh.edges()[iter.key()];

        meshTools::writeOBJ(pointStream, mesh.points()[e[0]]); vertI++;
        meshTools::writeOBJ(pointStream, mesh.points()[e[1]]); vertI++;
        pointStream<< "l " << vertI-1 << ' ' << vertI << nl;
    }
}



int main(int argc, char *argv[])
{
    argList::addNote
    (
        "for mesh debugging: write mesh as separate OBJ files"
    );

    timeSelector::addOptions();
    argList::addBoolOption
    (
        "patchFaces",
        "write patch faces edges"
    );
    argList::addBoolOption
    (
        "patchEdges",
        "write patch boundary edges"
    );
    argList::addOption
    (
        "cell",
        "int",
        "write points for the specified cell"
    );
    argList::addOption
    (
        "face",
        "int",
        "write specified face"
    );
    argList::addOption
    (
        "point",
        "int",
        "write specified point"
    );
    argList::addOption
    (
        "cellSet",
        "name",
        "write points for specified cellSet"
    );
    argList::addOption
    (
        "faceSet",
        "name",
        "write points for specified faceSet"
    );
    #include "addRegionOption.H"

    #include "setRootCase.H"
    #include "createTime.H"
    runTime.functionObjects().off();

    const bool patchFaces = args.optionFound("patchFaces");
    const bool patchEdges = args.optionFound("patchEdges");
    const bool doCell     = args.optionFound("cell");
    const bool doPoint    = args.optionFound("point");
    const bool doFace     = args.optionFound("face");
    const bool doCellSet  = args.optionFound("cellSet");
    const bool doFaceSet  = args.optionFound("faceSet");


    Info<< "Writing mesh objects as .obj files such that the object"
        << " numbering" << endl
        << "(for points, faces, cells) is consistent with"
        << " Foam numbering (starting from 0)." << endl << endl;

    instantList timeDirs = timeSelector::select0(runTime, args);

    #include "createNamedPolyMesh.H"

    forAll(timeDirs, timeI)
    {
        runTime.setTime(timeDirs[timeI], timeI);

        Info<< "Time = " << runTime.timeName() << endl;

        polyMesh::readUpdateState state = mesh.readUpdate();

        if (!timeI || state != polyMesh::UNCHANGED)
        {
            if (patchFaces)
            {
                writePatchFaces(mesh, runTime.timeName());
            }
            if (patchEdges)
            {
                writePatchBoundaryEdges(mesh, runTime.timeName());
            }
            if (doCell)
            {
                label celli = args.optionRead<label>("cell");

                writePoints(mesh, celli, runTime.timeName());
            }
            if (doPoint)
            {
                label pointi = args.optionRead<label>("point");

                writePointCells(mesh, pointi, runTime.timeName());
            }
            if (doFace)
            {
                label facei = args.optionRead<label>("face");

                fileName fName
                (
                    mesh.time().path()
                  / "meshPoints_"
                  + runTime.timeName()
                  + '_'
                  + name(facei)
                  + ".obj"
                );

                Info<< "Writing mesh points and edges to " << fName << endl;

                OFstream str(fName);

                const face& f = mesh.faces()[facei];

                meshTools::writeOBJ(str, faceList(1, f), mesh.points());
            }
            if (doCellSet)
            {
                const word setName = args["cellSet"];

                cellSet cells(mesh, setName);

                Info<< "Read " << cells.size() << " cells from set " << setName
                    << endl;

                writePoints(mesh, cells.toc(), runTime.timeName());
            }
            if (doFaceSet)
            {
                const word setName = args["faceSet"];

                faceSet faces(mesh, setName);

                Info<< "Read " << faces.size() << " faces from set " << setName
                    << endl;

                fileName fName
                (
                    mesh.time().path()
                  / "meshPoints_"
                  + runTime.timeName()
                  + '_'
                  + setName
                  + ".obj"
                );

                Info<< "Writing mesh points and edges to " << fName << endl;

                OFstream str(fName);

                meshTools::writeOBJ
                (
                    str,
                    mesh.faces(),
                    mesh.points(),
                    faces.toc()
                );
            }
            else if
            (
                !patchFaces
             && !patchEdges
             && !doCell
             && !doPoint
             && !doFace
             && !doCellSet
             && !doFaceSet
            )
            {
                // points & edges
                writePoints(mesh, runTime.timeName());

                // face centres
                writeFaceCentres(mesh, runTime.timeName());

                // cell centres
                writeCellCentres(mesh, runTime.timeName());

                // Patch face centres
                writePatchCentres(mesh, runTime.timeName());
            }
        }
        else
        {
            Info<< "No mesh." << endl;
        }

        Info<< nl << endl;
    }


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

    return 0;
}


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