writeMeshObj.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Application
25  writeMeshObj
26 
27 Description
28  For mesh debugging: writes mesh as three separate OBJ files which can
29  be viewed with e.g. javaview.
30 
31  meshPoints_XXX.obj : all points and edges as lines.
32  meshFaceCentres_XXX.obj : all face centres.
33  meshCellCentres_XXX.obj : all cell centres.
34 
35  patch_YYY_XXX.obj : all face centres of patch YYY
36 
37  Optional: - patch faces (as polygons) : patchFaces_YYY_XXX.obj
38  - non-manifold edges : patchEdges_YYY_XXX.obj
39 
40 \*---------------------------------------------------------------------------*/
41 
42 #include "argList.H"
43 #include "timeSelector.H"
44 #include "Time.H"
45 #include "polyMesh.H"
46 #include "OFstream.H"
47 #include "meshTools.H"
48 #include "cellSet.H"
49 #include "faceSet.H"
50 #include "SubField.H"
51 
52 using namespace Foam;
53 
54 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
55 
56 void writeOBJ(const point& pt, Ostream& os)
57 {
58  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
59 }
60 
61 // All edges of mesh
62 void writePoints(const polyMesh& mesh, const fileName& timeName)
63 {
64  fileName pointFile(mesh.time().path()/"meshPoints_" + timeName + ".obj");
65 
66  Info<< "Writing mesh points and edges to " << pointFile << endl;
67 
68  OFstream pointStream(pointFile);
69 
70  forAll(mesh.points(), pointi)
71  {
72  writeOBJ(mesh.points()[pointi], pointStream);
73  }
74 
75  forAll(mesh.edges(), edgeI)
76  {
77  const edge& e = mesh.edges()[edgeI];
78 
79  pointStream << "l " << e.start() + 1 << ' ' << e.end() + 1 << nl;
80  }
81 }
82 
83 
84 // Edges for subset of cells
85 void writePoints
86 (
87  const polyMesh& mesh,
88  const labelList& cellLabels,
89  const fileName& timeName
90 )
91 {
92  fileName fName(mesh.time().path()/"meshPoints_" + timeName + ".obj");
93 
94  Info<< "Writing mesh points and edges to " << fName << endl;
95 
96  OFstream str(fName);
97 
98  // OBJ file vertex
99  label vertI = 0;
100 
101  // From point to OBJ file vertex
102  Map<label> pointToObj(6*cellLabels.size());
103 
104  forAll(cellLabels, i)
105  {
106  label celli = cellLabels[i];
107 
108  const labelList& cEdges = mesh.cellEdges()[celli];
109 
110  forAll(cEdges, cEdgeI)
111  {
112  const edge& e = mesh.edges()[cEdges[cEdgeI]];
113 
114  label v0;
115 
116  Map<label>::iterator e0Fnd = pointToObj.find(e[0]);
117 
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  }
128 
129  label v1;
130 
131  Map<label>::iterator e1Fnd = pointToObj.find(e[1]);
132 
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  }
143 
144 
145  str << "l " << v0+1 << ' ' << v1+1 << nl;
146  }
147  }
148 }
149 
150 
151 // Edges of single cell
152 void writePoints
153 (
154  const polyMesh& mesh,
155  const label celli,
156  const fileName& timeName
157 )
158 {
159  fileName fName
160  (
161  mesh.time().path()
162  / "meshPoints_" + timeName + '_' + name(celli) + ".obj"
163  );
164 
165  Info<< "Writing mesh points and edges to " << fName << endl;
166 
167  OFstream pointStream(fName);
168 
169  const cell& cFaces = mesh.cells()[celli];
170 
171  meshTools::writeOBJ(pointStream, mesh.faces(), mesh.points(), cFaces);
172 }
173 
174 
175 
176 // All face centres
177 void writeFaceCentres(const polyMesh& mesh,const fileName& timeName)
178 {
179  fileName faceFile
180  (
181  mesh.time().path()
182  / "meshFaceCentres_" + timeName + ".obj"
183  );
184 
185  Info<< "Writing mesh face centres to " << faceFile << endl;
186 
187  OFstream faceStream(faceFile);
188 
189  forAll(mesh.faceCentres(), facei)
190  {
191  writeOBJ(mesh.faceCentres()[facei], faceStream);
192  }
193 }
194 
195 
196 void writeCellCentres(const polyMesh& mesh, const fileName& timeName)
197 {
198  fileName cellFile
199  (
200  mesh.time().path()/"meshCellCentres_" + timeName + ".obj"
201  );
202 
203  Info<< "Writing mesh cell centres to " << cellFile << endl;
204 
205  OFstream cellStream(cellFile);
206 
207  forAll(mesh.cellCentres(), celli)
208  {
209  writeOBJ(mesh.cellCentres()[celli], cellStream);
210  }
211 }
212 
213 
214 void writePatchCentres
215 (
216  const polyMesh& mesh,
217  const fileName& timeName
218 )
219 {
220  const polyBoundaryMesh& patches = mesh.boundaryMesh();
221 
223  {
224  const polyPatch& pp = patches[patchi];
225 
226  fileName faceFile
227  (
228  mesh.time().path()/"patch_" + pp.name() + '_' + timeName + ".obj"
229  );
230 
231  Info<< "Writing patch face centres to " << faceFile << endl;
232 
233  OFstream patchFaceStream(faceFile);
234 
235  forAll(pp.faceCentres(), facei)
236  {
237  writeOBJ(pp.faceCentres()[facei], patchFaceStream);
238  }
239  }
240 }
241 
242 
243 void writePatchFaces
244 (
245  const polyMesh& mesh,
246  const fileName& timeName
247 )
248 {
249  const polyBoundaryMesh& patches = mesh.boundaryMesh();
250 
252  {
253  const polyPatch& pp = patches[patchi];
254 
255  fileName faceFile
256  (
257  mesh.time().path()
258  / "patchFaces_" + pp.name() + '_' + timeName + ".obj"
259  );
260 
261  Info<< "Writing patch faces to " << faceFile << endl;
262 
263  OFstream patchFaceStream(faceFile);
264 
265  forAll(pp.localPoints(), pointi)
266  {
267  writeOBJ(pp.localPoints()[pointi], patchFaceStream);
268  }
269 
270  forAll(pp.localFaces(), facei)
271  {
272  const face& f = pp.localFaces()[facei];
273 
274  patchFaceStream<< 'f';
275 
276  forAll(f, fp)
277  {
278  patchFaceStream << ' ' << f[fp]+1;
279  }
280  patchFaceStream << nl;
281  }
282  }
283 }
284 
285 
286 void writePatchBoundaryEdges
287 (
288  const polyMesh& mesh,
289  const fileName& timeName
290 )
291 {
292  const polyBoundaryMesh& patches = mesh.boundaryMesh();
293 
295  {
296  const polyPatch& pp = patches[patchi];
297 
298  fileName edgeFile
299  (
300  mesh.time().path()
301  / "patchEdges_" + pp.name() + '_' + timeName + ".obj"
302  );
303 
304  Info<< "Writing patch edges to " << edgeFile << endl;
305 
306  OFstream patchEdgeStream(edgeFile);
307 
308  forAll(pp.localPoints(), pointi)
309  {
310  writeOBJ(pp.localPoints()[pointi], patchEdgeStream);
311  }
312 
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];
318 
319  patchEdgeStream<< "l " << e[0]+1 << ' ' << e[1]+1 << nl;
320  }
321  }
322  }
323 }
324 
325 
326 void writePointCells
327 (
328  const polyMesh& mesh,
329  const label pointi,
330  const fileName& timeName
331 )
332 {
333  const labelList& pCells = mesh.pointCells()[pointi];
334 
335  labelHashSet allEdges(6*pCells.size());
336 
337  forAll(pCells, i)
338  {
339  const labelList& cEdges = mesh.cellEdges()[pCells[i]];
340 
341  forAll(cEdges, i)
342  {
343  allEdges.insert(cEdges[i]);
344  }
345  }
346 
347 
348  fileName pFile
349  (
350  mesh.time().path()
351  / "pointEdges_" + timeName + '_' + name(pointi) + ".obj"
352  );
353 
354  Info<< "Writing pointEdges to " << pFile << endl;
355 
356  OFstream pointStream(pFile);
357 
358  label vertI = 0;
359 
360  forAllConstIter(labelHashSet, allEdges, iter)
361  {
362  const edge& e = mesh.edges()[iter.key()];
363 
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  }
368 }
369 
370 
371 
372 int main(int argc, char *argv[])
373 {
375  (
376  "for mesh debugging: write mesh as separate OBJ files"
377  );
378 
381  (
382  "patchFaces",
383  "write patch faces edges"
384  );
386  (
387  "patchEdges",
388  "write patch boundary edges"
389  );
391  (
392  "cell",
393  "int",
394  "write points for the specified cell"
395  );
397  (
398  "face",
399  "int",
400  "write specified face"
401  );
403  (
404  "point",
405  "int",
406  "write specified point"
407  );
409  (
410  "cellSet",
411  "name",
412  "write points for specified cellSet"
413  );
415  (
416  "faceSet",
417  "name",
418  "write points for specified faceSet"
419  );
420  #include "addRegionOption.H"
421 
422  #include "setRootCase.H"
424 
425  const bool patchFaces = args.optionFound("patchFaces");
426  const bool patchEdges = args.optionFound("patchEdges");
427  const bool doCell = args.optionFound("cell");
428  const bool doPoint = args.optionFound("point");
429  const bool doFace = args.optionFound("face");
430  const bool doCellSet = args.optionFound("cellSet");
431  const bool doFaceSet = args.optionFound("faceSet");
432 
433 
434  Info<< "Writing mesh objects as .obj files such that the object"
435  << " numbering" << endl
436  << "(for points, faces, cells) is consistent with"
437  << " Foam numbering (starting from 0)." << endl << endl;
438 
440 
441  #include "createNamedPolyMesh.H"
442 
443  forAll(timeDirs, timeI)
444  {
445  runTime.setTime(timeDirs[timeI], timeI);
446 
447  Info<< "Time = " << runTime.userTimeName() << endl;
448 
449  polyMesh::readUpdateState state = mesh.readUpdate();
450 
451  if (!timeI || state != polyMesh::UNCHANGED)
452  {
453  if (patchFaces)
454  {
455  writePatchFaces(mesh, runTime.name());
456  }
457  if (patchEdges)
458  {
459  writePatchBoundaryEdges(mesh, runTime.name());
460  }
461  if (doCell)
462  {
463  label celli = args.optionRead<label>("cell");
464 
465  writePoints(mesh, celli, runTime.name());
466  }
467  if (doPoint)
468  {
469  label pointi = args.optionRead<label>("point");
470 
471  writePointCells(mesh, pointi, runTime.name());
472  }
473  if (doFace)
474  {
475  label facei = args.optionRead<label>("face");
476 
477  fileName fName
478  (
479  mesh.time().path()
480  / "meshPoints_"
481  + runTime.name()
482  + '_'
483  + name(facei)
484  + ".obj"
485  );
486 
487  Info<< "Writing mesh points and edges to " << fName << endl;
488 
489  OFstream str(fName);
490 
491  const face& f = mesh.faces()[facei];
492 
493  meshTools::writeOBJ(str, faceList(1, f), mesh.points());
494  }
495  if (doCellSet)
496  {
497  const word setName = args["cellSet"];
498 
499  cellSet cells(mesh, setName);
500 
501  Info<< "Read " << cells.size() << " cells from set " << setName
502  << endl;
503 
504  writePoints(mesh, cells.toc(), runTime.name());
505  }
506  if (doFaceSet)
507  {
508  const word setName = args["faceSet"];
509 
510  faceSet faces(mesh, setName);
511 
512  Info<< "Read " << faces.size() << " faces from set " << setName
513  << endl;
514 
515  fileName fName
516  (
517  mesh.time().path()
518  / "meshPoints_"
519  + runTime.name()
520  + '_'
521  + setName
522  + ".obj"
523  );
524 
525  Info<< "Writing mesh points and edges to " << fName << endl;
526 
527  OFstream str(fName);
528 
530  (
531  str,
532  mesh.faces(),
533  mesh.points(),
534  faces.toc()
535  );
536  }
537  else if
538  (
539  !patchFaces
540  && !patchEdges
541  && !doCell
542  && !doPoint
543  && !doFace
544  && !doCellSet
545  && !doFaceSet
546  )
547  {
548  // points & edges
549  writePoints(mesh, runTime.name());
550 
551  // face centres
552  writeFaceCentres(mesh, runTime.name());
553 
554  // cell centres
555  writeCellCentres(mesh, runTime.name());
556 
557  // Patch face centres
558  writePatchCentres(mesh, runTime.name());
559  }
560  }
561  else
562  {
563  Info<< "No mesh." << endl;
564  }
565 
566  Info<< nl << endl;
567  }
568 
569 
570  Info<< "End\n" << endl;
571 
572  return 0;
573 }
574 
575 
576 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Output to file stream.
Definition: OFstream.H:86
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:57
label nEdges() const
Return number of edges in patch.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
label nInternalEdges() const
Number of internal edges.
const List< FaceType > & localFaces() const
Return patch faces addressing into local point list.
const labelListList & edgeFaces() const
Return edge-face addressing.
const Field< PointType > & localPoints() const
Return pointField of points in patch.
fileName path() const
Explicitly inherit path from TimePaths to disambiguate from.
Definition: TimePaths.H:138
const Cmpt & z() const
Definition: VectorI.H:87
const Cmpt & y() const
Definition: VectorI.H:81
const Cmpt & x() const
Definition: VectorI.H:75
static void addOption(const word &opt, const string &param="", const string &usage="")
Add to an option to validOptions with usage information.
Definition: argList.C:128
static void addNote(const string &)
Add extra notes for the usage information.
Definition: argList.C:159
T optionRead(const word &opt) const
Read a value from the named option.
Definition: argListI.H:193
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
A collection of cell labels.
Definition: cellSet.H:51
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:60
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:61
A list of face labels.
Definition: faceSet.H:51
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:76
A class for handling file names.
Definition: fileName.H:82
word name() const
Return file name (part beyond last /)
Definition: fileName.C:195
const Time & time() const
Return time.
const word & name() const
Return name.
Foam::polyBoundaryMesh.
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1373
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
readUpdateState readUpdate()
Update the mesh based on the mesh files saved in.
Definition: polyMeshIO.C:99
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:405
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1360
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:70
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:276
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const labelListList & cellEdges() const
const vectorField & faceCentres() const
const labelListList & pointCells() const
const vectorField & cellCentres() const
const cellList & cells() const
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
A class for handling words, derived from string.
Definition: word.H:62
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
static instantList timeDirs
Definition: globalFoam.H:44
const cellShapeList & cells
const fvPatchList & patches
word timeName
Definition: getTimeIndex.H:3
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of point.
Definition: meshTools.C:203
Namespace for OpenFOAM.
const doubleScalar e
Definition: doubleScalar.H:105
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
messageStream Info
List< face > faceList
Definition: faceListFwd.H:43
static const char nl
Definition: Ostream.H:260
labelList f(nPoints)
Foam::argList args(argc, argv)