initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / mesh / conversion / writeMeshObj / writeMeshObj.C
blobe1d571e0816f4f146decba57af55579348a03a2f
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Description
26     For mesh debugging: writes mesh as three separate OBJ files which can
27     be viewed with e.g. javaview.
29     meshPoints_XXX.obj : all points and edges as lines.
30     meshFaceCentres_XXX.obj : all face centres.
31     meshCellCentres_XXX.obj : all cell centres.
33     patch_YYY_XXX.obj : all face centres of patch YYY
35     Optional: patch faces (as polygons) : patchFaces_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"
48 using namespace Foam;
50 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
52 void writeOBJ(const point& pt, Ostream& os)
54     os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << endl;
57 // All edges of mesh
58 void writePoints(const polyMesh& mesh, const fileName& timeName)
60     label vertI = 0;
62     fileName pointFile(mesh.time().path()/"meshPoints_" + timeName + ".obj");
64     Info << "Writing mesh points and edges to " << pointFile << endl;
66     OFstream pointStream(pointFile);
68     forAll(mesh.points(), pointI)
69     {
70         writeOBJ(mesh.points()[pointI], pointStream);
71         vertI++;
72     }
74     forAll(mesh.edges(), edgeI)
75     {
76         const edge& e = mesh.edges()[edgeI];
78         pointStream << "l " << e.start() + 1 << ' ' << e.end() + 1
79             << endl;
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 << endl;
281         }
282     }
286 void writePointCells
288     const polyMesh& mesh,
289     const label pointI,
290     const fileName& timeName
293     const labelList& pCells = mesh.pointCells()[pointI];
295     labelHashSet allEdges(6*pCells.size());
297     forAll(pCells, i)
298     {
299         const labelList& cEdges = mesh.cellEdges()[pCells[i]];
301         forAll(cEdges, i)
302         {
303             allEdges.insert(cEdges[i]);
304         }
305     }
308     fileName pFile
309     (
310         mesh.time().path()
311       / "pointEdges_" + timeName + '_' + name(pointI) + ".obj"
312     );
314     Info << "Writing pointEdges to " << pFile << endl;
316     OFstream pointStream(pFile);
318     label vertI = 0;
320     for
321     (
322         labelHashSet::const_iterator iter = allEdges.begin();
323         iter != allEdges.end();
324         ++iter
325     )
326     {
327         const edge& e = mesh.edges()[iter.key()];
329         meshTools::writeOBJ(pointStream, mesh.points()[e[0]]); vertI++;
330         meshTools::writeOBJ(pointStream, mesh.points()[e[1]]); vertI++;
331         pointStream<< "l " << vertI-1 << ' ' << vertI << nl;
332     }
336 // Main program:
338 int main(int argc, char *argv[])
340     timeSelector::addOptions();
341     argList::validOptions.insert("patchFaces", "");
342     argList::validOptions.insert("cell", "cellI");
343     argList::validOptions.insert("face", "faceI");
344     argList::validOptions.insert("point", "pointI");
345     argList::validOptions.insert("cellSet", "setName");
346     argList::validOptions.insert("faceSet", "setName");
347 #   include "addRegionOption.H"
349 #   include "setRootCase.H"
350 #   include "createTime.H"
351     runTime.functionObjects().off();
353     bool patchFaces = args.optionFound("patchFaces");
354     bool doCell     = args.optionFound("cell");
355     bool doPoint    = args.optionFound("point");
356     bool doFace     = args.optionFound("face");
357     bool doCellSet  = args.optionFound("cellSet");
358     bool doFaceSet  = args.optionFound("faceSet");
361     Info<< "Writing mesh objects as .obj files such that the object"
362         << " numbering" << endl
363         << "(for points, faces, cells) is consistent with"
364         << " Foam numbering (starting from 0)." << endl << endl;
366     instantList timeDirs = timeSelector::select0(runTime, args);
368 #   include "createNamedPolyMesh.H"
370     forAll(timeDirs, timeI)
371     {
372         runTime.setTime(timeDirs[timeI], timeI);
374         Info<< "Time = " << runTime.timeName() << endl;
376         polyMesh::readUpdateState state = mesh.readUpdate();
378         if (!timeI || state != polyMesh::UNCHANGED)
379         {
380             if (patchFaces)
381             {
382                 writePatchFaces(mesh, runTime.timeName());
383             }
384             else if (doCell)
385             {
386                 label cellI = args.optionRead<label>("cell");
388                 writePoints(mesh, cellI, runTime.timeName());
389             }
390             else if (doPoint)
391             {
392                 label pointI = args.optionRead<label>("point");
394                 writePointCells(mesh, pointI, runTime.timeName());
395             }
396             else if (doFace)
397             {
398                 label faceI = args.optionRead<label>("face");
400                 fileName fName
401                 (
402                     mesh.time().path()
403                   / "meshPoints_"
404                   + runTime.timeName()
405                   + '_'
406                   + name(faceI)
407                   + ".obj"
408                 );
410                 Info<< "Writing mesh points and edges to " << fName << endl;
412                 OFstream str(fName);
414                 const face& f = mesh.faces()[faceI];
416                 meshTools::writeOBJ(str, faceList(1, f), mesh.points());
417             }
418             else if (doCellSet)
419             {
420                 word setName(args.option("cellSet"));
422                 cellSet cells(mesh, setName);
424                 Info<< "Read " << cells.size() << " cells from set " << setName
425                     << endl;
427                 writePoints(mesh, cells.toc(), runTime.timeName());
429             }
430             else if (doFaceSet)
431             {
432                 word setName(args.option("faceSet"));
434                 faceSet faces(mesh, setName);
436                 Info<< "Read " << faces.size() << " faces from set " << setName
437                     << endl;
439                 fileName fName
440                 (
441                     mesh.time().path()
442                   / "meshPoints_"
443                   + runTime.timeName()
444                   + '_'
445                   + setName
446                   + ".obj"
447                 );
449                 Info << "Writing mesh points and edges to " << fName << endl;
451                 OFstream str(fName);
453                 meshTools::writeOBJ
454                 (
455                     str,
456                     mesh.faces(),
457                     mesh.points(),
458                     faces.toc()
459                 );
460             }
461             else
462             {
463                 // points & edges
464                 writePoints(mesh, runTime.timeName());
466                 // face centres
467                 writeFaceCentres(mesh, runTime.timeName());
469                 // cell centres
470                 writeCellCentres(mesh, runTime.timeName());
472                 // Patch face centres
473                 writePatchCentres(mesh, runTime.timeName());
474             }
475         }
476         else
477         {
478             Info << "No mesh." << endl;
479         }
481         Info << nl << endl;
482     }
485     Info << "End\n" << endl;
487     return 0;
491 // ************************************************************************* //