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  label vertI = 0;
71 
72  fileName pointFile(mesh.time().path()/"meshPoints_" + timeName + ".obj");
73 
74  Info<< "Writing mesh points and edges to " << pointFile << endl;
75 
76  OFstream pointStream(pointFile);
77 
78  forAll(mesh.points(), pointi)
79  {
80  writeOBJ(mesh.points()[pointi], pointStream);
81  vertI++;
82  }
83 
84  forAll(mesh.edges(), edgeI)
85  {
86  const edge& e = mesh.edges()[edgeI];
87 
88  pointStream << "l " << e.start() + 1 << ' ' << e.end() + 1 << nl;
89  }
90 }
91 
92 
93 // Edges for subset of cells
94 void writePoints
95 (
96  const polyMesh& mesh,
97  const labelList& cellLabels,
98  const fileName& timeName
99 )
100 {
101  fileName fName(mesh.time().path()/"meshPoints_" + timeName + ".obj");
102 
103  Info<< "Writing mesh points and edges to " << fName << endl;
104 
105  OFstream str(fName);
106 
107  // OBJ file vertex
108  label vertI = 0;
109 
110  // From point to OBJ file vertex
111  Map<label> pointToObj(6*cellLabels.size());
112 
113  forAll(cellLabels, i)
114  {
115  label celli = cellLabels[i];
116 
117  const labelList& cEdges = mesh.cellEdges()[celli];
118 
119  forAll(cEdges, cEdgeI)
120  {
121  const edge& e = mesh.edges()[cEdges[cEdgeI]];
122 
123  label v0;
124 
125  const auto e0Fnd = pointToObj.cfind(e[0]);
126 
127  if (e0Fnd.found())
128  {
129  v0 = *e0Fnd;
130  }
131  else
132  {
133  v0 = vertI++;
134  meshTools::writeOBJ(str, mesh.points()[e[0]]);
135  pointToObj.insert(e[0], v0);
136  }
137 
138  label v1;
139 
140  const auto e1Fnd = pointToObj.cfind(e[1]);
141 
142  if (e1Fnd.found())
143  {
144  v1 = *e1Fnd;
145  }
146  else
147  {
148  v1 = vertI++;
149  meshTools::writeOBJ(str, mesh.points()[e[1]]);
150  pointToObj.insert(e[1], v1);
151  }
152 
153 
154  str << "l " << v0+1 << ' ' << v1+1 << nl;
155  }
156  }
157 }
158 
159 
160 // Edges of single cell
161 void writePoints
162 (
163  const polyMesh& mesh,
164  const label celli,
165  const fileName& timeName
166 )
167 {
168  fileName fName
169  (
170  mesh.time().path()
171  / "meshPoints_" + timeName + '_' + name(celli) + ".obj"
172  );
173 
174  Info<< "Writing mesh points and edges to " << fName << endl;
175 
176  OFstream pointStream(fName);
177 
178  const cell& cFaces = mesh.cells()[celli];
179 
180  meshTools::writeOBJ(pointStream, mesh.faces(), mesh.points(), cFaces);
181 }
182 
183 
184 
185 // All face centres
186 void writeFaceCentres(const polyMesh& mesh,const fileName& timeName)
187 {
188  fileName faceFile
189  (
190  mesh.time().path()
191  / "meshFaceCentres_" + timeName + ".obj"
192  );
193 
194  Info<< "Writing mesh face centres to " << faceFile << endl;
195 
196  OFstream faceStream(faceFile);
197 
198  forAll(mesh.faceCentres(), facei)
199  {
200  writeOBJ(mesh.faceCentres()[facei], faceStream);
201  }
202 }
203 
204 
205 void writeCellCentres(const polyMesh& mesh, const fileName& timeName)
206 {
207  fileName cellFile
208  (
209  mesh.time().path()/"meshCellCentres_" + timeName + ".obj"
210  );
211 
212  Info<< "Writing mesh cell centres to " << cellFile << endl;
213 
214  OFstream cellStream(cellFile);
215 
216  forAll(mesh.cellCentres(), celli)
217  {
218  writeOBJ(mesh.cellCentres()[celli], cellStream);
219  }
220 }
221 
222 
223 void writePatchCentres
224 (
225  const polyMesh& mesh,
226  const fileName& timeName
227 )
228 {
230 
231  forAll(patches, patchi)
232  {
233  const polyPatch& pp = patches[patchi];
234 
235  fileName faceFile
236  (
237  mesh.time().path()/"patch_" + pp.name() + '_' + timeName + ".obj"
238  );
239 
240  Info<< "Writing patch face centres to " << faceFile << endl;
241 
242  OFstream patchFaceStream(faceFile);
243 
244  forAll(pp.faceCentres(), facei)
245  {
246  writeOBJ(pp.faceCentres()[facei], patchFaceStream);
247  }
248  }
249 }
250 
251 
252 void writePatchFaces
253 (
254  const polyMesh& mesh,
255  const fileName& timeName
256 )
257 {
259 
260  forAll(patches, patchi)
261  {
262  const polyPatch& pp = patches[patchi];
263 
264  fileName faceFile
265  (
266  mesh.time().path()
267  / "patchFaces_" + pp.name() + '_' + timeName + ".obj"
268  );
269 
270  Info<< "Writing patch faces to " << faceFile << endl;
271 
272  OFstream patchFaceStream(faceFile);
273 
274  forAll(pp.localPoints(), pointi)
275  {
276  writeOBJ(pp.localPoints()[pointi], patchFaceStream);
277  }
278 
279  forAll(pp.localFaces(), facei)
280  {
281  const face& f = pp.localFaces()[facei];
282 
283  patchFaceStream<< 'f';
284 
285  forAll(f, fp)
286  {
287  patchFaceStream << ' ' << f[fp]+1;
288  }
289  patchFaceStream << nl;
290  }
291  }
292 }
293 
294 
295 void writePatchBoundaryEdges
296 (
297  const polyMesh& mesh,
298  const fileName& timeName
299 )
300 {
302 
303  forAll(patches, patchi)
304  {
305  const polyPatch& pp = patches[patchi];
306 
307  fileName edgeFile
308  (
309  mesh.time().path()
310  / "patchEdges_" + pp.name() + '_' + timeName + ".obj"
311  );
312 
313  Info<< "Writing patch edges to " << edgeFile << endl;
314 
315  OFstream patchEdgeStream(edgeFile);
316 
317  forAll(pp.localPoints(), pointi)
318  {
319  writeOBJ(pp.localPoints()[pointi], patchEdgeStream);
320  }
321 
322  for (label edgeI = pp.nInternalEdges(); edgeI < pp.nEdges(); edgeI++)
323  {
324  if (pp.edgeFaces()[edgeI].size() == 1)
325  {
326  const edge& e = pp.edges()[edgeI];
327 
328  patchEdgeStream<< "l " << e[0]+1 << ' ' << e[1]+1 << nl;
329  }
330  }
331  }
332 }
333 
334 
335 void writePointCells
336 (
337  const polyMesh& mesh,
338  const label pointi,
339  const fileName& timeName
340 )
341 {
342  const labelList& pCells = mesh.pointCells()[pointi];
343 
344  labelHashSet allEdges(6*pCells.size());
345 
346  forAll(pCells, i)
347  {
348  const labelList& cEdges = mesh.cellEdges()[pCells[i]];
349 
350  allEdges.insert(cEdges);
351  }
352 
353 
354  fileName pFile
355  (
356  mesh.time().path()
357  / "pointEdges_" + timeName + '_' + name(pointi) + ".obj"
358  );
359 
360  Info<< "Writing pointEdges to " << pFile << endl;
361 
362  OFstream pointStream(pFile);
363 
364  label vertI = 0;
365 
366  for (const label edgei : allEdges)
367  {
368  const edge& e = mesh.edges()[edgei];
369 
370  meshTools::writeOBJ(pointStream, mesh.points()[e[0]]); vertI++;
371  meshTools::writeOBJ(pointStream, mesh.points()[e[1]]); vertI++;
372  pointStream<< "l " << vertI-1 << ' ' << vertI << nl;
373  }
374 }
375 
376 
377 
378 int main(int argc, char *argv[])
379 {
381  (
382  "For mesh debugging: write mesh as separate OBJ files"
383  );
384 
387  (
388  "patchFaces",
389  "Write patch faces edges"
390  );
392  (
393  "patchEdges",
394  "Write patch boundary edges"
395  );
397  (
398  "cell",
399  "cellId",
400  "Write points for the specified cell"
401  );
403  (
404  "face",
405  "faceId",
406  "Write specified face"
407  );
409  (
410  "point",
411  "pointId",
412  "Write specified point"
413  );
415  (
416  "cellSet",
417  "name",
418  "Write points for specified cellSet"
419  );
421  (
422  "faceSet",
423  "name",
424  "Write points for specified faceSet"
425  );
426  #include "addRegionOption.H"
427 
428  argList::noFunctionObjects(); // Never use function objects
429 
430  #include "setRootCase.H"
431  #include "createTime.H"
432 
433  const bool patchFaces = args.found("patchFaces");
434  const bool patchEdges = args.found("patchEdges");
435  const bool doCell = args.found("cell");
436  const bool doPoint = args.found("point");
437  const bool doFace = args.found("face");
438  const bool doCellSet = args.found("cellSet");
439  const bool doFaceSet = args.found("faceSet");
440 
441 
442  Info<< "Writing mesh objects as .obj files such that the object"
443  << " numbering" << endl
444  << "(for points, faces, cells) is consistent with"
445  << " Foam numbering (starting from 0)." << endl << endl;
446 
448 
449  #include "createNamedPolyMesh.H"
450 
451  forAll(timeDirs, timeI)
452  {
453  runTime.setTime(timeDirs[timeI], timeI);
454 
455  Info<< "Time = " << runTime.timeName() << endl;
456 
458 
459  if (!timeI || state != polyMesh::UNCHANGED)
460  {
461  if (patchFaces)
462  {
463  writePatchFaces(mesh, runTime.timeName());
464  }
465  if (patchEdges)
466  {
467  writePatchBoundaryEdges(mesh, runTime.timeName());
468  }
469  if (doCell)
470  {
471  const label celli = args.get<label>("cell");
472 
473  writePoints(mesh, celli, runTime.timeName());
474  }
475  if (doPoint)
476  {
477  const label pointi = args.get<label>("point");
478 
479  writePointCells(mesh, pointi, runTime.timeName());
480  }
481  if (doFace)
482  {
483  const label facei = args.get<label>("face");
484 
485  fileName fName
486  (
487  mesh.time().path()
488  / "meshPoints_"
489  + runTime.timeName()
490  + '_'
491  + name(facei)
492  + ".obj"
493  );
494 
495  Info<< "Writing mesh points and edges to " << fName << endl;
496 
497  OFstream str(fName);
498 
499  const face& f = mesh.faces()[facei];
500 
502  }
503  if (doCellSet)
504  {
505  const word setName = args["cellSet"];
506 
507  cellSet cells(mesh, setName);
508 
509  Info<< "Read " << cells.size() << " cells from set " << setName
510  << endl;
511 
512  writePoints(mesh, cells.toc(), runTime.timeName());
513  }
514  if (doFaceSet)
515  {
516  const word setName = args["faceSet"];
517 
518  faceSet faces(mesh, setName);
519 
520  Info<< "Read " << faces.size() << " faces from set " << setName
521  << endl;
522 
523  fileName fName
524  (
525  mesh.time().path()
526  / "meshPoints_"
527  + runTime.timeName()
528  + '_'
529  + setName
530  + ".obj"
531  );
532 
533  Info<< "Writing mesh points and edges to " << fName << endl;
534 
535  OFstream str(fName);
536 
538  (
539  str,
540  mesh.faces(),
541  mesh.points(),
542  faces.toc()
543  );
544  }
545  else if
546  (
547  !patchFaces
548  && !patchEdges
549  && !doCell
550  && !doPoint
551  && !doFace
552  && !doCellSet
553  && !doFaceSet
554  )
555  {
556  // points & edges
557  writePoints(mesh, runTime.timeName());
558 
559  // face centres
560  writeFaceCentres(mesh, runTime.timeName());
561 
562  // cell centres
563  writeCellCentres(mesh, runTime.timeName());
564 
565  // Patch face centres
566  writePatchCentres(mesh, runTime.timeName());
567  }
568  }
569  else
570  {
571  Info<< "No mesh." << endl;
572  }
573 
574  Info<< nl << endl;
575  }
576 
577 
578  Info<< "End\n" << endl;
579 
580  return 0;
581 }
582 
583 
584 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::polyMesh::points
virtual const pointField & points() const
Return raw points.
Definition: polyMesh.C:1069
meshTools.H
Foam::Vector::x
const Cmpt & x() const
Access to the vector x component.
Definition: VectorI.H:73
Foam::PrimitivePatch::edgeFaces
const labelListList & edgeFaces() const
Return edge-face addressing.
Definition: PrimitivePatch.C:262
Foam::pointFile
Inserts points at locations specified in a pointFile into the surfaces to be conformed to of the conf...
Definition: pointFile.H:55
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fileName
A class for handling file names.
Definition: fileName.H:73
SubField.H
Foam::polyBoundaryMesh
A polyBoundaryMesh is a polyPatch list with additional search methods and registered IO.
Definition: polyBoundaryMesh.H:63
Foam::PrimitivePatch::edges
const edgeList & edges() const
Return list of edges, address into LOCAL point list.
Definition: PrimitivePatch.C:183
Foam::meshTools::writeOBJ
void writeOBJ(Ostream &os, const point &pt)
Write obj representation of a point.
Definition: meshTools.C:203
Foam::PrimitivePatch::nEdges
label nEdges() const
Number of edges in patch.
Definition: PrimitivePatch.H:322
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
Foam::primitiveMesh::cells
const cellList & cells() const
Definition: primitiveMeshCells.C:138
Foam::edge
An edge is a list of two point labels. The functionality it provides supports the discretisation on a...
Definition: edge.H:63
Foam::argList::addNote
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:412
Foam::Map< label >
Foam::faceSet
A list of face labels.
Definition: faceSet.H:51
Foam::polyMesh::boundaryMesh
const polyBoundaryMesh & boundaryMesh() const
Return boundary mesh.
Definition: polyMesh.H:444
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::primitiveMesh::edges
const edgeList & edges() const
Return mesh edges. Uses calcEdges.
Definition: primitiveMeshEdges.C:505
polyMesh.H
Foam::HashSet< label, Hash< label > >
Foam::argList::get
T get(const label index) const
Get a value from the argument at index.
Definition: argListI.H:278
Foam::Vector::z
const Cmpt & z() const
Access to the vector z component.
Definition: VectorI.H:85
createNamedPolyMesh.H
Required Variables.
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:77
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
OFstream.H
Foam::fvMesh::readUpdate
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:648
Foam::argList::noFunctionObjects
static void noFunctionObjects(bool addWithOption=false)
Remove '-noFunctionObjects' option and ignore any occurrences.
Definition: argList.C:473
Foam::primitiveMesh::cellEdges
const labelListList & cellEdges() const
Definition: primitiveMeshCellEdges.C:115
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::polyPatch
A patch is a list of labels that address the faces in the global face list.
Definition: polyPatch.H:68
argList.H
faceSet.H
addRegionOption.H
Foam::polyMesh::UNCHANGED
Definition: polyMesh.H:92
timeName
word timeName
Definition: getTimeIndex.H:3
os
OBJstream os(runTime.globalPath()/outputName)
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::PrimitivePatch::localFaces
const List< face_type > & localFaces() const
Return patch faces addressing into local point list.
Definition: PrimitivePatch.C:317
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::cellSet
A collection of cell labels.
Definition: cellSet.H:51
Foam::PrimitivePatch::localPoints
const Field< point_type > & localPoints() const
Return pointField of points in patch.
Definition: PrimitivePatch.C:359
Foam::OFstream
Output to file stream, using an OSstream.
Definition: OFstream.H:53
Foam::polyMesh::readUpdateState
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:90
Foam::PrimitivePatch::nInternalEdges
label nInternalEdges() const
Number of internal edges.
Definition: PrimitivePatch.C:214
Foam::argList::addBoolOption
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:324
Time.H
setRootCase.H
Foam::Vector::y
const Cmpt & y() const
Access to the vector y component.
Definition: VectorI.H:79
Foam::primitiveMesh::cellCentres
const vectorField & cellCentres() const
Definition: primitiveMeshCellCentresAndVols.C:84
Foam::polyMesh::faces
virtual const faceList & faces() const
Return raw faces.
Definition: polyMesh.C:1094
Foam::nl
constexpr char nl
Definition: Ostream.H:404
Foam::Time::path
fileName path() const
Return path.
Definition: Time.H:358
Foam::polyPatch::faceCentres
const vectorField::subField faceCentres() const
Return face centres.
Definition: polyPatch.C:321
f
labelList f(nPoints)
Foam::faceList
List< face > faceList
A List of faces.
Definition: faceListFwd.H:47
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Add timeSelector options to argList::validOptions.
Definition: timeSelector.C:102
Foam::Vector< scalar >
Foam::List< label >
Foam::Time::setTime
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:1003
Foam::primitiveMesh::faceCentres
const vectorField & faceCentres() const
Definition: primitiveMeshFaceCentresAndAreas.C:77
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
Foam::primitiveMesh::pointCells
const labelListList & pointCells() const
Definition: primitiveMeshPointCells.C:110
timeSelector.H
createTime.H
patches
const polyBoundaryMesh & patches
Definition: convertProcessorPatches.H:65
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::fvMesh::time
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:280
Foam::face
A face is a list of labels corresponding to mesh vertices.
Definition: face.H:72
Foam::patchIdentifier::name
const word & name() const noexcept
The patch name.
Definition: patchIdentifier.H:135
cells
const cellShapeList & cells
Definition: gmvOutputHeader.H:3
cellSet.H
Foam::Ostream
An Ostream is an abstract base class for all output systems (streams, files, token lists,...
Definition: Ostream.H:56
Foam::timeSelector::select0
static instantList select0(Time &runTime, const argList &args)
Definition: timeSelector.C:235
Foam::argList::addOption
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:335
args
Foam::argList args(argc, argv)
Foam::cell
A cell is defined as a list of faces with extra functionality.
Definition: cell.H:54
Foam::argList::found
bool found(const word &optName) const
Return true if the named option is found.
Definition: argListI.H:178