Transferred copyright to the OpenFOAM Foundation
[OpenFOAM-2.0.x.git] / applications / utilities / mesh / conversion / writeMeshObj / writeMeshObj.C
blob068e8415abe5fa2f04a8683ef5fcdb69019636be
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
11     OpenFOAM is free software: you can redistribute it and/or modify it
12     under the terms of the GNU General Public License as published by
13     the Free Software Foundation, either version 3 of the License, or
14     (at your option) any later version.
16     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
24 Description
25     For mesh debugging: writes mesh as three separate OBJ files which can
26     be viewed with e.g. javaview.
28     meshPoints_XXX.obj : all points and edges as lines.
29     meshFaceCentres_XXX.obj : all face centres.
30     meshCellCentres_XXX.obj : all cell centres.
32     patch_YYY_XXX.obj : all face centres of patch YYY
34     Optional: - patch faces (as polygons) : patchFaces_YYY_XXX.obj
35               - non-manifold edges : patchEdges_YYY_XXX.obj
37 \*---------------------------------------------------------------------------*/
39 #include "argList.H"
40 #include "timeSelector.H"
41 #include "Time.H"
42 #include "polyMesh.H"
43 #include "OFstream.H"
44 #include "meshTools.H"
45 #include "cellSet.H"
46 #include "faceSet.H"
47 #include "SubField.H"
49 using namespace Foam;
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 void writeOBJ(const point& pt, Ostream& os)
55     os  << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
58 // All edges of mesh
59 void writePoints(const polyMesh& mesh, const fileName& timeName)
61     label vertI = 0;
63     fileName pointFile(mesh.time().path()/"meshPoints_" + timeName + ".obj");
65     Info<< "Writing mesh points and edges to " << pointFile << endl;
67     OFstream pointStream(pointFile);
69     forAll(mesh.points(), pointI)
70     {
71         writeOBJ(mesh.points()[pointI], pointStream);
72         vertI++;
73     }
75     forAll(mesh.edges(), edgeI)
76     {
77         const edge& e = mesh.edges()[edgeI];
79         pointStream << "l " << e.start() + 1 << ' ' << e.end() + 1 << nl;
80     }
84 // Edges for subset of cells
85 void writePoints
87     const polyMesh& mesh,
88     const labelList& cellLabels,
89     const fileName& timeName
92     fileName fName(mesh.time().path()/"meshPoints_" + timeName + ".obj");
94     Info<< "Writing mesh points and edges to " << fName << endl;
96     OFstream str(fName);
98     // OBJ file vertex
99     label vertI = 0;
101     // From point to OBJ file vertex
102     Map<label> pointToObj(6*cellLabels.size());
104     forAll(cellLabels, i)
105     {
106         label cellI = cellLabels[i];
108         const labelList& cEdges = mesh.cellEdges()[cellI];
110         forAll(cEdges, cEdgeI)
111         {
112             const edge& e = mesh.edges()[cEdges[cEdgeI]];
114             label v0;
116             Map<label>::iterator e0Fnd = pointToObj.find(e[0]);
118             if (e0Fnd == pointToObj.end())
119             {
120                 meshTools::writeOBJ(str, mesh.points()[e[0]]);
121                 v0 = vertI++;
122                 pointToObj.insert(e[0], v0);
123             }
124             else
125             {
126                 v0 = e0Fnd();
127             }
129             label v1;
131             Map<label>::iterator e1Fnd = pointToObj.find(e[1]);
133             if (e1Fnd == pointToObj.end())
134             {
135                 meshTools::writeOBJ(str, mesh.points()[e[1]]);
136                 v1 = vertI++;
137                 pointToObj.insert(e[1], v1);
138             }
139             else
140             {
141                 v1 = e1Fnd();
142             }
145             str << "l " << v0+1 << ' ' << v1+1 << nl;
146         }
147     }
151 // Edges of single cell
152 void writePoints
154     const polyMesh& mesh,
155     const label cellI,
156     const fileName& timeName
159     fileName fName
160     (
161         mesh.time().path()
162       / "meshPoints_" + timeName + '_' + name(cellI) + ".obj"
163     );
165     Info<< "Writing mesh points and edges to " << fName << endl;
167     OFstream pointStream(fName);
169     const cell& cFaces = mesh.cells()[cellI];
171     meshTools::writeOBJ(pointStream, mesh.faces(), mesh.points(), cFaces);
176 // All face centres
177 void writeFaceCentres(const polyMesh& mesh,const fileName& timeName)
179     fileName faceFile
180     (
181         mesh.time().path()
182       / "meshFaceCentres_" + timeName + ".obj"
183     );
185     Info<< "Writing mesh face centres to " << faceFile << endl;
187     OFstream faceStream(faceFile);
189     forAll(mesh.faceCentres(), faceI)
190     {
191         writeOBJ(mesh.faceCentres()[faceI], faceStream);
192     }
196 void writeCellCentres(const polyMesh& mesh, const fileName& timeName)
198     fileName cellFile
199     (
200         mesh.time().path()/"meshCellCentres_" + timeName + ".obj"
201     );
203     Info<< "Writing mesh cell centres to " << cellFile << endl;
205     OFstream cellStream(cellFile);
207     forAll(mesh.cellCentres(), cellI)
208     {
209         writeOBJ(mesh.cellCentres()[cellI], cellStream);
210     }
214 void writePatchCentres
216     const polyMesh& mesh,
217     const fileName& timeName
220     const polyBoundaryMesh& patches = mesh.boundaryMesh();
222     forAll(patches, patchI)
223     {
224         const polyPatch& pp = patches[patchI];
226         fileName faceFile
227         (
228             mesh.time().path()/"patch_" + pp.name() + '_' + timeName + ".obj"
229         );
231         Info<< "Writing patch face centres to " << faceFile << endl;
233         OFstream patchFaceStream(faceFile);
235         forAll(pp.faceCentres(), faceI)
236         {
237             writeOBJ(pp.faceCentres()[faceI], patchFaceStream);
238         }
239     }
243 void writePatchFaces
245     const polyMesh& mesh,
246     const fileName& timeName
249     const polyBoundaryMesh& patches = mesh.boundaryMesh();
251     forAll(patches, patchI)
252     {
253         const polyPatch& pp = patches[patchI];
255         fileName faceFile
256         (
257             mesh.time().path()
258           / "patchFaces_" + pp.name() + '_' + timeName + ".obj"
259         );
261         Info<< "Writing patch faces to " << faceFile << endl;
263         OFstream patchFaceStream(faceFile);
265         forAll(pp.localPoints(), pointI)
266         {
267             writeOBJ(pp.localPoints()[pointI], patchFaceStream);
268         }
270         forAll(pp.localFaces(), faceI)
271         {
272             const face& f = pp.localFaces()[faceI];
274             patchFaceStream<< 'f';
276             forAll(f, fp)
277             {
278                 patchFaceStream << ' ' << f[fp]+1;
279             }
280             patchFaceStream << nl;
281         }
282     }
286 void writePatchBoundaryEdges
288     const polyMesh& mesh,
289     const fileName& timeName
292     const polyBoundaryMesh& patches = mesh.boundaryMesh();
294     forAll(patches, patchI)
295     {
296         const polyPatch& pp = patches[patchI];
298         fileName edgeFile
299         (
300             mesh.time().path()
301           / "patchEdges_" + pp.name() + '_' + timeName + ".obj"
302         );
304         Info<< "Writing patch edges to " << edgeFile << endl;
306         OFstream patchEdgeStream(edgeFile);
308         forAll(pp.localPoints(), pointI)
309         {
310             writeOBJ(pp.localPoints()[pointI], patchEdgeStream);
311         }
313         for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
314         {
315             if (pp.edgeFaces()[edgeI].size() == 1)
316             {
317                 const edge& e = pp.edges()[edgeI];
319                 patchEdgeStream<< "l " << e[0]+1 << ' ' << e[1]+1 << nl;
320             }
321         }
322     }
326 void writePointCells
328     const polyMesh& mesh,
329     const label pointI,
330     const fileName& timeName
333     const labelList& pCells = mesh.pointCells()[pointI];
335     labelHashSet allEdges(6*pCells.size());
337     forAll(pCells, i)
338     {
339         const labelList& cEdges = mesh.cellEdges()[pCells[i]];
341         forAll(cEdges, i)
342         {
343             allEdges.insert(cEdges[i]);
344         }
345     }
348     fileName pFile
349     (
350         mesh.time().path()
351       / "pointEdges_" + timeName + '_' + name(pointI) + ".obj"
352     );
354     Info<< "Writing pointEdges to " << pFile << endl;
356     OFstream pointStream(pFile);
358     label vertI = 0;
360     forAllConstIter(labelHashSet, allEdges, iter)
361     {
362         const edge& e = mesh.edges()[iter.key()];
364         meshTools::writeOBJ(pointStream, mesh.points()[e[0]]); vertI++;
365         meshTools::writeOBJ(pointStream, mesh.points()[e[1]]); vertI++;
366         pointStream<< "l " << vertI-1 << ' ' << vertI << nl;
367     }
371 // Main program:
373 int main(int argc, char *argv[])
375     argList::addNote
376     (
377         "for mesh debugging: write mesh as separate OBJ files"
378     );
380     timeSelector::addOptions();
381     argList::addBoolOption
382     (
383         "patchFaces",
384         "write patch faces edges"
385     );
386     argList::addBoolOption
387     (
388         "patchEdges",
389         "write patch boundary edges"
390     );
391     argList::addOption
392     (
393         "cell",
394         "int",
395         "write points for the specified cell"
396     );
397     argList::addOption
398     (
399         "face",
400         "int",
401         "write specified face"
402     );
403     argList::addOption
404     (
405         "point",
406         "int",
407         "write specified point"
408     );
409     argList::addOption
410     (
411         "cellSet",
412         "name",
413         "write points for specified cellSet"
414     );
415     argList::addOption
416     (
417         "faceSet",
418         "name",
419         "write points for specified faceSet"
420     );
421     #include "addRegionOption.H"
423     #include "setRootCase.H"
424     #include "createTime.H"
425     runTime.functionObjects().off();
427     const bool patchFaces = args.optionFound("patchFaces");
428     const bool patchEdges = args.optionFound("patchEdges");
429     const bool doCell     = args.optionFound("cell");
430     const bool doPoint    = args.optionFound("point");
431     const bool doFace     = args.optionFound("face");
432     const bool doCellSet  = args.optionFound("cellSet");
433     const bool doFaceSet  = args.optionFound("faceSet");
436     Info<< "Writing mesh objects as .obj files such that the object"
437         << " numbering" << endl
438         << "(for points, faces, cells) is consistent with"
439         << " Foam numbering (starting from 0)." << endl << endl;
441     instantList timeDirs = timeSelector::select0(runTime, args);
443     #include "createNamedPolyMesh.H"
445     forAll(timeDirs, timeI)
446     {
447         runTime.setTime(timeDirs[timeI], timeI);
449         Info<< "Time = " << runTime.timeName() << endl;
451         polyMesh::readUpdateState state = mesh.readUpdate();
453         if (!timeI || state != polyMesh::UNCHANGED)
454         {
455             if (patchFaces)
456             {
457                 writePatchFaces(mesh, runTime.timeName());
458             }
459             if (patchEdges)
460             {
461                 writePatchBoundaryEdges(mesh, runTime.timeName());
462             }
463             if (doCell)
464             {
465                 label cellI = args.optionRead<label>("cell");
467                 writePoints(mesh, cellI, runTime.timeName());
468             }
469             if (doPoint)
470             {
471                 label pointI = args.optionRead<label>("point");
473                 writePointCells(mesh, pointI, runTime.timeName());
474             }
475             if (doFace)
476             {
477                 label faceI = args.optionRead<label>("face");
479                 fileName fName
480                 (
481                     mesh.time().path()
482                   / "meshPoints_"
483                   + runTime.timeName()
484                   + '_'
485                   + name(faceI)
486                   + ".obj"
487                 );
489                 Info<< "Writing mesh points and edges to " << fName << endl;
491                 OFstream str(fName);
493                 const face& f = mesh.faces()[faceI];
495                 meshTools::writeOBJ(str, faceList(1, f), mesh.points());
496             }
497             if (doCellSet)
498             {
499                 const word setName = args["cellSet"];
501                 cellSet cells(mesh, setName);
503                 Info<< "Read " << cells.size() << " cells from set " << setName
504                     << endl;
506                 writePoints(mesh, cells.toc(), runTime.timeName());
507             }
508             if (doFaceSet)
509             {
510                 const word setName = args["faceSet"];
512                 faceSet faces(mesh, setName);
514                 Info<< "Read " << faces.size() << " faces from set " << setName
515                     << endl;
517                 fileName fName
518                 (
519                     mesh.time().path()
520                   / "meshPoints_"
521                   + runTime.timeName()
522                   + '_'
523                   + setName
524                   + ".obj"
525                 );
527                 Info<< "Writing mesh points and edges to " << fName << endl;
529                 OFstream str(fName);
531                 meshTools::writeOBJ
532                 (
533                     str,
534                     mesh.faces(),
535                     mesh.points(),
536                     faces.toc()
537                 );
538             }
539             else if
540             (
541                 !patchFaces
542              && !patchEdges
543              && !doCell
544              && !doPoint
545              && !doFace
546              && !doCellSet
547              && !doFaceSet
548             )
549             {
550                 // points & edges
551                 writePoints(mesh, runTime.timeName());
553                 // face centres
554                 writeFaceCentres(mesh, runTime.timeName());
556                 // cell centres
557                 writeCellCentres(mesh, runTime.timeName());
559                 // Patch face centres
560                 writePatchCentres(mesh, runTime.timeName());
561             }
562         }
563         else
564         {
565             Info<< "No mesh." << endl;
566         }
568         Info<< nl << endl;
569     }
572     Info<< "End\n" << endl;
574     return 0;
578 // ************************************************************************* //