1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
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
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/>.
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 \*---------------------------------------------------------------------------*/
40 #include "timeSelector.H"
44 #include "meshTools.H"
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 void writeOBJ(const point& pt, Ostream& os)
55 os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
59 void writePoints(const polyMesh& mesh, const fileName& timeName)
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)
71 writeOBJ(mesh.points()[pointI], pointStream);
75 forAll(mesh.edges(), edgeI)
77 const edge& e = mesh.edges()[edgeI];
79 pointStream << "l " << e.start() + 1 << ' ' << e.end() + 1 << nl;
84 // Edges for subset of cells
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;
101 // From point to OBJ file vertex
102 Map<label> pointToObj(6*cellLabels.size());
104 forAll(cellLabels, i)
106 label cellI = cellLabels[i];
108 const labelList& cEdges = mesh.cellEdges()[cellI];
110 forAll(cEdges, cEdgeI)
112 const edge& e = mesh.edges()[cEdges[cEdgeI]];
116 Map<label>::iterator e0Fnd = pointToObj.find(e[0]);
118 if (e0Fnd == pointToObj.end())
120 meshTools::writeOBJ(str, mesh.points()[e[0]]);
122 pointToObj.insert(e[0], v0);
131 Map<label>::iterator e1Fnd = pointToObj.find(e[1]);
133 if (e1Fnd == pointToObj.end())
135 meshTools::writeOBJ(str, mesh.points()[e[1]]);
137 pointToObj.insert(e[1], v1);
145 str << "l " << v0+1 << ' ' << v1+1 << nl;
151 // Edges of single cell
154 const polyMesh& mesh,
156 const fileName& timeName
162 / "meshPoints_" + timeName + '_' + name(cellI) + ".obj"
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);
177 void writeFaceCentres(const polyMesh& mesh,const fileName& timeName)
182 / "meshFaceCentres_" + timeName + ".obj"
185 Info<< "Writing mesh face centres to " << faceFile << endl;
187 OFstream faceStream(faceFile);
189 forAll(mesh.faceCentres(), faceI)
191 writeOBJ(mesh.faceCentres()[faceI], faceStream);
196 void writeCellCentres(const polyMesh& mesh, const fileName& timeName)
200 mesh.time().path()/"meshCellCentres_" + timeName + ".obj"
203 Info<< "Writing mesh cell centres to " << cellFile << endl;
205 OFstream cellStream(cellFile);
207 forAll(mesh.cellCentres(), cellI)
209 writeOBJ(mesh.cellCentres()[cellI], cellStream);
214 void writePatchCentres
216 const polyMesh& mesh,
217 const fileName& timeName
220 const polyBoundaryMesh& patches = mesh.boundaryMesh();
222 forAll(patches, patchI)
224 const polyPatch& pp = patches[patchI];
228 mesh.time().path()/"patch_" + pp.name() + '_' + timeName + ".obj"
231 Info<< "Writing patch face centres to " << faceFile << endl;
233 OFstream patchFaceStream(faceFile);
235 forAll(pp.faceCentres(), faceI)
237 writeOBJ(pp.faceCentres()[faceI], patchFaceStream);
245 const polyMesh& mesh,
246 const fileName& timeName
249 const polyBoundaryMesh& patches = mesh.boundaryMesh();
251 forAll(patches, patchI)
253 const polyPatch& pp = patches[patchI];
258 / "patchFaces_" + pp.name() + '_' + timeName + ".obj"
261 Info<< "Writing patch faces to " << faceFile << endl;
263 OFstream patchFaceStream(faceFile);
265 forAll(pp.localPoints(), pointI)
267 writeOBJ(pp.localPoints()[pointI], patchFaceStream);
270 forAll(pp.localFaces(), faceI)
272 const face& f = pp.localFaces()[faceI];
274 patchFaceStream<< 'f';
278 patchFaceStream << ' ' << f[fp]+1;
280 patchFaceStream << nl;
286 void writePatchBoundaryEdges
288 const polyMesh& mesh,
289 const fileName& timeName
292 const polyBoundaryMesh& patches = mesh.boundaryMesh();
294 forAll(patches, patchI)
296 const polyPatch& pp = patches[patchI];
301 / "patchEdges_" + pp.name() + '_' + timeName + ".obj"
304 Info<< "Writing patch edges to " << edgeFile << endl;
306 OFstream patchEdgeStream(edgeFile);
308 forAll(pp.localPoints(), pointI)
310 writeOBJ(pp.localPoints()[pointI], patchEdgeStream);
313 for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
315 if (pp.edgeFaces()[edgeI].size() == 1)
317 const edge& e = pp.edges()[edgeI];
319 patchEdgeStream<< "l " << e[0]+1 << ' ' << e[1]+1 << nl;
328 const polyMesh& mesh,
330 const fileName& timeName
333 const labelList& pCells = mesh.pointCells()[pointI];
335 labelHashSet allEdges(6*pCells.size());
339 const labelList& cEdges = mesh.cellEdges()[pCells[i]];
343 allEdges.insert(cEdges[i]);
351 / "pointEdges_" + timeName + '_' + name(pointI) + ".obj"
354 Info<< "Writing pointEdges to " << pFile << endl;
356 OFstream pointStream(pFile);
360 forAllConstIter(labelHashSet, allEdges, iter)
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;
373 int main(int argc, char *argv[])
377 "for mesh debugging: write mesh as separate OBJ files"
380 timeSelector::addOptions();
381 argList::addBoolOption
384 "write patch faces edges"
386 argList::addBoolOption
389 "write patch boundary edges"
395 "write points for the specified cell"
401 "write specified face"
407 "write specified point"
413 "write points for specified cellSet"
419 "write points for specified faceSet"
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)
447 runTime.setTime(timeDirs[timeI], timeI);
449 Info<< "Time = " << runTime.timeName() << endl;
451 polyMesh::readUpdateState state = mesh.readUpdate();
453 if (!timeI || state != polyMesh::UNCHANGED)
457 writePatchFaces(mesh, runTime.timeName());
461 writePatchBoundaryEdges(mesh, runTime.timeName());
465 label cellI = args.optionRead<label>("cell");
467 writePoints(mesh, cellI, runTime.timeName());
471 label pointI = args.optionRead<label>("point");
473 writePointCells(mesh, pointI, runTime.timeName());
477 label faceI = args.optionRead<label>("face");
489 Info<< "Writing mesh points and edges to " << fName << endl;
493 const face& f = mesh.faces()[faceI];
495 meshTools::writeOBJ(str, faceList(1, f), mesh.points());
499 const word setName = args["cellSet"];
501 cellSet cells(mesh, setName);
503 Info<< "Read " << cells.size() << " cells from set " << setName
506 writePoints(mesh, cells.toc(), runTime.timeName());
510 const word setName = args["faceSet"];
512 faceSet faces(mesh, setName);
514 Info<< "Read " << faces.size() << " faces from set " << setName
527 Info<< "Writing mesh points and edges to " << fName << endl;
551 writePoints(mesh, runTime.timeName());
554 writeFaceCentres(mesh, runTime.timeName());
557 writeCellCentres(mesh, runTime.timeName());
559 // Patch face centres
560 writePatchCentres(mesh, runTime.timeName());
565 Info<< "No mesh." << endl;
572 Info<< "End\n" << endl;
578 // ************************************************************************* //