writeMeshObj.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2018 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 Application
28  writeMeshObj
29 
30 Group
31  grpMeshConversionUtilities
32 
33 Description
34  For mesh debugging: writes mesh as three separate OBJ files.
35 
36  meshPoints_XXX.obj : all points and edges as lines.
37  meshFaceCentres_XXX.obj : all face centres.
38  meshCellCentres_XXX.obj : all cell centres.
39 
40  patch_YYY_XXX.obj : all face centres of patch YYY
41 
42  Optional: - patch faces (as polygons) : patchFaces_YYY_XXX.obj
43  - non-manifold edges : patchEdges_YYY_XXX.obj
44 
45 \*---------------------------------------------------------------------------*/
46 
47 #include "argList.H"
48 #include "timeSelector.H"
49 #include "Time.H"
50 #include "polyMesh.H"
51 #include "OFstream.H"
52 #include "meshTools.H"
53 #include "cellSet.H"
54 #include "faceSet.H"
55 #include "SubField.H"
56 
57 using namespace Foam;
58 
59 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
60 
61 void writeOBJ(const point& pt, Ostream& os)
62 {
63  os << "v " << pt.x() << ' ' << pt.y() << ' ' << pt.z() << nl;
64 }
65 
66 
67 // All edges of mesh
68 void writePoints(const polyMesh& mesh, const fileName& timeName)
69 {
70  fileName pointFile(mesh.time().path()/"meshPoints_" + timeName + ".obj");
71 
72  Info<< "Writing mesh points and edges to " << pointFile << endl;
73 
74  OFstream pointStream(pointFile);
75 
76  forAll(mesh.points(), pointi)
77  {
78  writeOBJ(mesh.points()[pointi], pointStream);
79  }
80 
81  forAll(mesh.edges(), edgeI)
82  {
83  const edge& e = mesh.edges()[edgeI];
84 
85  pointStream << "l " << e.start() + 1 << ' ' << e.end() + 1 << nl;
86  }
87 }
88 
89 
90 // Edges for subset of cells
91 void writePoints
92 (
93  const polyMesh& mesh,
94  const labelList& cellLabels,
95  const fileName& timeName
96 )
97 {
98  fileName fName(mesh.time().path()/"meshPoints_" + timeName + ".obj");
99 
100  Info<< "Writing mesh points and edges to " << fName << endl;
101 
102  OFstream str(fName);
103 
104  // OBJ file vertex
105  label vertI = 0;
106 
107  // From point to OBJ file vertex
108  Map<label> pointToObj(6*cellLabels.size());
109 
110  forAll(cellLabels, i)
111  {
112  label celli = cellLabels[i];
113 
114  const labelList& cEdges = mesh.cellEdges()[celli];
115 
116  forAll(cEdges, cEdgeI)
117  {
118  const edge& e = mesh.edges()[cEdges[cEdgeI]];
119 
120  label v0;
121 
122  const auto e0Fnd = pointToObj.cfind(e[0]);
123 
124  if (e0Fnd.good())
125  {
126  v0 = *e0Fnd;
127  }
128  else
129  {
130  v0 = vertI++;
131  meshTools::writeOBJ(str, mesh.points()[e[0]]);
132  pointToObj.insert(e[0], v0);
133  }
134 
135  label v1;
136 
137  const auto e1Fnd = pointToObj.cfind(e[1]);
138 
139  if (e1Fnd.good())
140  {
141  v1 = *e1Fnd;
142  }
143  else
144  {
145  v1 = vertI++;
146  meshTools::writeOBJ(str, mesh.points()[e[1]]);
147  pointToObj.insert(e[1], v1);
148  }
149 
150 
151  str << "l " << v0+1 << ' ' << v1+1 << nl;
152  }
153  }
154 }
155 
156 
157 // Edges of single cell
158 void writePoints
159 (
160  const polyMesh& mesh,
161  const label celli,
162  const fileName& timeName
163 )
164 {
165  fileName fName
166  (
167  mesh.time().path()
168  / "meshPoints_" + timeName + '_' + name(celli) + ".obj"
169  );
170 
171  Info<< "Writing mesh points and edges to " << fName << endl;
172 
173  OFstream pointStream(fName);
174 
175  const cell& cFaces = mesh.cells()[celli];
176 
177  meshTools::writeOBJ(pointStream, mesh.faces(), mesh.points(), cFaces);
178 }
179 
180 
181 
182 // All face centres
183 void writeFaceCentres(const polyMesh& mesh,const fileName& timeName)
184 {
185  fileName faceFile
186  (
187  mesh.time().path()
188  / "meshFaceCentres_" + timeName + ".obj"
189  );
190 
191  Info<< "Writing mesh face centres to " << faceFile << endl;
192 
193  OFstream faceStream(faceFile);
194 
195  forAll(mesh.faceCentres(), facei)
196  {
197  writeOBJ(mesh.faceCentres()[facei], faceStream);
198  }
199 }
200 
201 
202 void writeCellCentres(const polyMesh& mesh, const fileName& timeName)
203 {
204  fileName cellFile
205  (
206  mesh.time().path()/"meshCellCentres_" + timeName + ".obj"
207  );
208 
209  Info<< "Writing mesh cell centres to " << cellFile << endl;
210 
211  OFstream cellStream(cellFile);
212 
213  forAll(mesh.cellCentres(), celli)
214  {
215  writeOBJ(mesh.cellCentres()[celli], cellStream);
216  }
217 }
218 
219 
220 void writePatchCentres
221 (
222  const polyMesh& mesh,
223  const fileName& timeName
224 )
225 {
227 
228  forAll(patches, patchi)
229  {
230  const polyPatch& pp = patches[patchi];
231 
232  fileName faceFile
233  (
234  mesh.time().path()/"patch_" + pp.name() + '_' + timeName + ".obj"
235  );
236 
237  Info<< "Writing patch face centres to " << faceFile << endl;
238 
239  OFstream patchFaceStream(faceFile);
240 
241  forAll(pp.faceCentres(), facei)
242  {
243  writeOBJ(pp.faceCentres()[facei], patchFaceStream);
244  }
245  }
246 }
247 
248 
249 void writePatchFaces
250 (
251  const polyMesh& mesh,
252  const fileName& timeName
253 )
254 {
256 
257  forAll(patches, patchi)
258  {
259  const polyPatch& pp = patches[patchi];
260 
261  fileName faceFile
262  (
263  mesh.time().path()
264  / "patchFaces_" + pp.name() + '_' + timeName + ".obj"
265  );
266 
267  Info<< "Writing patch faces to " << faceFile << endl;
268 
269  OFstream patchFaceStream(faceFile);
270 
271  forAll(pp.localPoints(), pointi)
272  {
273  writeOBJ(pp.localPoints()[pointi], patchFaceStream);
274  }
275 
276  forAll(pp.localFaces(), facei)
277  {
278  const face& f = pp.localFaces()[facei];
279 
280  patchFaceStream<< 'f';
281 
282  forAll(f, fp)
283  {
284  patchFaceStream << ' ' << f[fp]+1;
285  }
286  patchFaceStream << nl;
287  }
288  }
289 }
290 
291 
292 void writePatchBoundaryEdges
293 (
294  const polyMesh& mesh,
295  const fileName& timeName
296 )
297 {
299 
300  forAll(patches, patchi)
301  {
302  const polyPatch& pp = patches[patchi];
303 
304  fileName edgeFile
305  (
306  mesh.time().path()
307  / "patchEdges_" + pp.name() + '_' + timeName + ".obj"
308  );
309 
310  Info<< "Writing patch edges to " << edgeFile << endl;
311 
312  OFstream patchEdgeStream(edgeFile);
313 
314  forAll(pp.localPoints(), pointi)
315  {
316  writeOBJ(pp.localPoints()[pointi], patchEdgeStream);
317  }
318 
319  for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
320  {
321  if (pp.edgeFaces()[edgeI].size() == 1)
322  {
323  const edge& e = pp.edges()[edgeI];
324 
325  patchEdgeStream<< "l " << e[0]+1 << ' ' << e[1]+1 << nl;
326  }
327  }
328  }
329 }
330 
331 
332 void writePointCells
333 (
334  const polyMesh& mesh,
335  const label pointi,
336  const fileName& timeName
337 )
338 {
339  const labelList& pCells = mesh.pointCells()[pointi];
340 
341  labelHashSet allEdges(6*pCells.size());
342 
343  forAll(pCells, i)
344  {
345  const labelList& cEdges = mesh.cellEdges()[pCells[i]];
346 
347  allEdges.insert(cEdges);
348  }
349 
350 
351  fileName pFile
352  (
353  mesh.time().path()
354  / "pointEdges_" + timeName + '_' + name(pointi) + ".obj"
355  );
356 
357  Info<< "Writing pointEdges to " << pFile << endl;
358 
359  OFstream pointStream(pFile);
360 
361  label vertI = 0;
362 
363  for (const label edgei : allEdges)
364  {
365  const edge& e = mesh.edges()[edgei];
366 
367  meshTools::writeOBJ(pointStream, mesh.points()[e[0]]); vertI++;
368  meshTools::writeOBJ(pointStream, mesh.points()[e[1]]); vertI++;
369  pointStream<< "l " << vertI-1 << ' ' << vertI << nl;
370  }
371 }
372 
373 
374 
375 int main(int argc, char *argv[])
376 {
378  (
379  "For mesh debugging: write mesh as separate OBJ files"
380  );
381 
384  (
385  "patchFaces",
386  "Write patch faces edges"
387  );
389  (
390  "patchEdges",
391  "Write patch boundary edges"
392  );
394  (
395  "cell",
396  "cellId",
397  "Write points for the specified cell"
398  );
400  (
401  "face",
402  "faceId",
403  "Write specified face"
404  );
406  (
407  "point",
408  "pointId",
409  "Write specified point"
410  );
412  (
413  "cellSet",
414  "name",
415  "Write points for specified cellSet"
416  );
418  (
419  "faceSet",
420  "name",
421  "Write points for specified faceSet"
422  );
423  #include "addRegionOption.H"
424 
425  argList::noFunctionObjects(); // Never use function objects
426 
427  #include "setRootCase.H"
428  #include "createTime.H"
429 
430  const bool patchFaces = args.found("patchFaces");
431  const bool patchEdges = args.found("patchEdges");
432  const bool doCell = args.found("cell");
433  const bool doPoint = args.found("point");
434  const bool doFace = args.found("face");
435  const bool doCellSet = args.found("cellSet");
436  const bool doFaceSet = args.found("faceSet");
437 
438 
439  Info<< "Writing mesh objects as .obj files such that the object"
440  << " numbering" << endl
441  << "(for points, faces, cells) is consistent with"
442  << " Foam numbering (starting from 0)." << endl << endl;
443 
445 
446  #include "createNamedPolyMesh.H"
447 
448  forAll(timeDirs, timeI)
449  {
450  runTime.setTime(timeDirs[timeI], timeI);
451 
452  Info<< "Time = " << runTime.timeName() << endl;
453 
455 
456  if (!timeI || state != polyMesh::UNCHANGED)
457  {
458  if (patchFaces)
459  {
460  writePatchFaces(mesh, runTime.timeName());
461  }
462  if (patchEdges)
463  {
464  writePatchBoundaryEdges(mesh, runTime.timeName());
465  }
466  if (doCell)
467  {
468  const label celli = args.get<label>("cell");
469 
470  writePoints(mesh, celli, runTime.timeName());
471  }
472  if (doPoint)
473  {
474  const label pointi = args.get<label>("point");
475 
476  writePointCells(mesh, pointi, runTime.timeName());
477  }
478  if (doFace)
479  {
480  const label facei = args.get<label>("face");
481 
482  fileName fName
483  (
484  mesh.time().path()
485  / "meshPoints_"
486  + runTime.timeName()
487  + '_'
488  + name(facei)
489  + ".obj"
490  );
491 
492  Info<< "Writing mesh points and edges to " << fName << endl;
493 
494  OFstream str(fName);
495 
496  const face& f = mesh.faces()[facei];
497 
499  }
500  if (doCellSet)
501  {
502  const word setName = args["cellSet"];
503 
504  cellSet cells(mesh, setName);
505 
506  Info<< "Read " << cells.size() << " cells from set " << setName
507  << endl;
508 
509  writePoints(mesh, cells.toc(), runTime.timeName());
510  }
511  if (doFaceSet)
512  {
513  const word setName = args["faceSet"];
514 
515  faceSet faces(mesh, setName);
516 
517  Info<< "Read " << faces.size() << " faces from set " << setName
518  << endl;
519 
520  fileName fName
521  (
522  mesh.time().path()
523  / "meshPoints_"
524  + runTime.timeName()
525  + '_'
526  + setName
527  + ".obj"
528  );
529 
530  Info<< "Writing mesh points and edges to " << fName << endl;
531 
532  OFstream str(fName);
533 
535  (
536  str,
537  mesh.faces(),
538  mesh.points(),
539  faces.toc()
540  );
541  }
542  else if
543  (
544  !patchFaces
545  && !patchEdges
546  && !doCell
547  && !doPoint
548  && !doFace
549  && !doCellSet
550  && !doFaceSet
551  )
552  {
553  // points & edges
554  writePoints(mesh, runTime.timeName());
555 
556  // face centres
557  writeFaceCentres(mesh, runTime.timeName());
558 
559  // cell centres
560  writeCellCentres(mesh, runTime.timeName());
561 
562  // Patch face centres
563  writePatchCentres(mesh, runTime.timeName());
564  }
565  }
566  else
567  {
568  Info<< "No mesh." << endl;
569  }
570 
571  Info<< nl << endl;
572  }
573 
574 
575  Info<< "End\n" << endl;
576 
577  return 0;
578 }
579 
580 
581 // ************************************************************************* //
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
const labelListList & cellEdges() const
static void noFunctionObjects(bool addWithOption=false)
Remove &#39;-noFunctionObjects&#39; option and ignore any occurrences.
Definition: argList.C:561
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:476
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
fileName path() const
Return path = rootPath/caseName. Same as TimePaths::path()
Definition: Time.H:503
A class for handling file names.
Definition: fileName.H:72
A list of face labels.
Definition: faceSet.H:47
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:68
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Output to file stream as an OSstream, normally using std::ofstream for the actual output...
Definition: OFstream.H:71
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:50
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:529
label nInternalEdges() const
Number of internal edges.
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:388
const cellList & cells() const
Required Classes.
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:196
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:360
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1063
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:286
List< face > faceList
List of faces.
Definition: faceListFwd.H:39
word timeName
Definition: getTimeIndex.H:3
Required Classes.
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:670
const cellShapeList & cells
An edge is a list of two vertex labels. This can correspond to a directed graph edge or an edge on a ...
Definition: edge.H:59
const polyBoundaryMesh & boundaryMesh() const noexcept
Return boundary mesh.
Definition: polyMesh.H:609
A class for handling words, derived from Foam::string.
Definition: word.H:63
const Field< point_type > & faceCentres() const
Return face centres for patch.
static void addOption(const word &optName, const string &param="", const string &usage="", bool advanced=false)
Add an option to validOptions with usage information.
Definition: argList.C:399
const labelListList & edgeFaces() const
Return edge-face addressing.
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:902
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1088
const vectorField & cellCentres() const
const labelListList & pointCells() const
A polyBoundaryMesh is a polyPatch list with registered IO, a reference to the associated polyMesh...
An Ostream is an abstract base class for all output systems (streams, files, token lists...
Definition: Ostream.H:56
label nEdges() const
Number of edges in patch.
A Vector of values with scalar precision, where scalar is float/double depending on the compilation f...
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
OBJstream os(runTime.globalPath()/outputName)
labelList f(nPoints)
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options and also set the runTime to the first i...
Definition: timeSelector.C:260
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
const vectorField & faceCentres() const
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:271
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:53
A collection of cell labels.
Definition: cellSet.H:47
const polyBoundaryMesh & patches
messageStream Info
Information stream (stdout output on master, null elsewhere)
Inserts points at locations specified in a pointFile into the surfaces to be conformed to of the conf...
Definition: pointFile.H:51
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:75
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:69
Foam::argList args(argc, argv)
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:171
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:101
uindirectPrimitivePatch pp(UIndirectList< face >(mesh.faces(), faceLabels), mesh.points())
Namespace for OpenFOAM.