foamToEnsight.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-2022 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  foamToEnsight
26 
27 Description
28  Translates OpenFOAM data to EnSight format.
29 
30  An Ensight part is created for the internalMesh and for each patch.
31 
32 Usage
33  \b foamToEnsight [OPTION]
34 
35  Options:
36  - \par -ascii
37  Write Ensight data in ASCII format instead of "C Binary"
38 
39  - \par -patches patchList
40  Specify particular patches to write.
41  Specifying an empty list suppresses writing the internalMesh.
42 
43  - \par -noPatches
44  Suppress writing any patches.
45 
46  - \par -faceZones zoneList
47  Specify faceZones to write, with wildcards
48 
49  - \par -cellZone zoneName
50  Specify single cellZone to write (not lagrangian)
51 
52  Note:
53  Parallel support for cloud data is not supported
54  - writes to \a EnSight directory to avoid collisions with
55  foamToEnsightParts
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #include "argList.H"
60 #include "timeSelector.H"
61 #include "IOobjectList.H"
62 #include "IOmanip.H"
63 #include "OFstream.H"
64 
65 #include "volFields.H"
66 
67 #include "labelIOField.H"
68 #include "scalarIOField.H"
69 #include "tensorIOField.H"
70 
71 #include "ensightMesh.H"
72 #include "ensightField.H"
73 
75 #include "ensightCloudField.H"
76 
77 #include "fvc.H"
78 
79 #include "cellSet.H"
80 #include "fvMeshSubset.H"
81 
82 using namespace Foam;
83 
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86 bool inFileNameList
87 (
88  const fileNameList& nameList,
89  const word& name
90 )
91 {
92  forAll(nameList, i)
93  {
94  if (nameList[i] == name)
95  {
96  return true;
97  }
98  }
99 
100  return false;
101 }
102 
103 
104 
105 int main(int argc, char *argv[])
106 {
108  #include "addRegionOption.H"
109 
111  (
112  "ascii",
113  "write in ASCII format instead of 'C Binary'"
114  );
116  (
117  "nodeValues",
118  "write values in nodes"
119  );
121  (
122  "noPatches",
123  "suppress writing any patches"
124  );
126  (
127  "patches",
128  "wordReList",
129  "specify particular patches to write - eg '(outlet \"inlet.*\")'. "
130  "An empty list suppresses writing the internalMesh."
131  );
133  (
134  "faceZones",
135  "wordReList",
136  "specify faceZones to write - eg '( slice \"mfp-.*\" )'."
137  );
139  (
140  "fields",
141  "wordReList",
142  "specify fields to export (all by default) - eg '( \"U.*\" )'."
143  );
145  (
146  "cellZone",
147  "word",
148  "specify cellZone to write"
149  );
150 
151  #include "setRootCase.H"
152 
153  // Check options
154  const bool binary = !args.optionFound("ascii");
155  const bool nodeValues = args.optionFound("nodeValues");
156 
157  #include "createTime.H"
158 
159  instantList Times = timeSelector::select0(runTime, args);
160 
161  #include "createNamedMesh.H"
162 
163  // Mesh instance (region0 gets filtered out)
164  fileName regionPrefix = "";
165 
167  {
168  regionPrefix = regionName;
169  }
170 
171  const label nVolFieldTypes = 5;
172  const word volFieldTypes[] =
173  {
179  };
180 
181  // Path to EnSight directory at case level only
182  // - For parallel cases, data only written from master
183  fileName ensightDir = args.rootPath()/args.globalCaseName()/"EnSight";
184 
185  if (Pstream::master())
186  {
187  if (isDir(ensightDir))
188  {
189  rmDir(ensightDir);
190  }
191 
192  mkDir(ensightDir);
193  }
194 
195  // Start of case file header output
196  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
197 
198  const word prepend = args.globalCaseName() + '.';
199 
200  OFstream *ensightCaseFilePtr = nullptr;
201  if (Pstream::master())
202  {
203  fileName caseFileName = prepend + "case";
204  Info<< nl << "write case: " << caseFileName.c_str() << endl;
205 
206  // the case file is always ASCII
207  ensightCaseFilePtr = new OFstream
208  (
209  ensightDir/caseFileName,
211  );
212 
213  *ensightCaseFilePtr
214  << "FORMAT" << nl
215  << "type: ensight gold" << nl << nl;
216  }
217 
218  OFstream& ensightCaseFile = *ensightCaseFilePtr;
219 
220  // Construct the EnSight mesh
221  const bool selectedPatches = args.optionFound("patches");
222  wordReList patchPatterns;
223  if (selectedPatches)
224  {
225  patchPatterns = wordReList(args.optionLookup("patches")());
226  }
227  const bool selectedZones = args.optionFound("faceZones");
228  wordReList zonePatterns;
229  if (selectedZones)
230  {
231  zonePatterns = wordReList(args.optionLookup("faceZones")());
232  }
233 
234  const bool selectedFields = args.optionFound("fields");
235  wordReList fieldPatterns;
236  if (selectedFields)
237  {
238  fieldPatterns = wordReList(args.optionLookup("fields")());
239  }
240 
241  word cellZoneName;
242  const bool doCellZone = args.optionReadIfPresent("cellZone", cellZoneName);
243 
244  fvMeshSubset meshSubsetter(mesh);
245  if (doCellZone)
246  {
247  Info<< "Converting cellZone " << cellZoneName
248  << " only (puts outside faces into patch "
249  << mesh.boundaryMesh()[0].name()
250  << ")" << endl;
251  const cellZone& cz = mesh.cellZones()[cellZoneName];
252  cellSet c0(mesh, "c0", labelHashSet(cz));
253  meshSubsetter.setLargeCellSubset(c0, 0);
254  }
255 
256  ensightMesh eMesh
257  (
258  (
259  meshSubsetter.hasSubMesh()
260  ? meshSubsetter.subMesh()
261  : meshSubsetter.baseMesh()
262  ),
263  args.optionFound("noPatches"),
264  selectedPatches,
265  patchPatterns,
266  selectedZones,
267  zonePatterns,
268  binary
269  );
270 
271  // Set Time to the last time before looking for the lagrangian objects
272  runTime.setTime(Times.last(), Times.size()-1);
273 
274  IOobjectList objects(mesh, runTime.name());
275 
276  #include "checkMeshMoving.H"
277 
278  if (meshMoving)
279  {
280  Info<< "Detected a moving mesh (multiple polyMesh/points files)."
281  << " Writing meshes for every timestep." << endl;
282  }
283 
284 
285  wordHashSet allCloudNames;
286  if (Pstream::master())
287  {
288  word geomFileName = prepend + "0000";
289 
290  // test pre check variable if there is a moving mesh
291  if (meshMoving)
292  {
293  geomFileName = prepend + "****";
294  }
295 
296  ensightCaseFile
297  << "GEOMETRY" << nl
298  << "model: 1 "
299  << (geomFileName + ".mesh").c_str() << nl;
300  }
301 
302  // Identify if lagrangian data exists at each time, and add clouds
303  // to the 'allCloudNames' hash set
304  forAll(Times, timeI)
305  {
306  runTime.setTime(Times[timeI], timeI);
307 
308  fileNameList cloudDirs = readDir
309  (
310  runTime.timePath()/regionPrefix/cloud::prefix,
312  );
313 
314  forAll(cloudDirs, cloudI)
315  {
316  IOobjectList cloudObjs
317  (
318  mesh,
319  runTime.name(),
320  cloud::prefix/cloudDirs[cloudI]
321  );
322 
323  IOobject* positionsPtr = cloudObjs.lookup(word("positions"));
324 
325  if (positionsPtr)
326  {
327  allCloudNames.insert(cloudDirs[cloudI]);
328  }
329  }
330  }
331 
332  HashTable<HashTable<word>> allCloudFields;
333  forAllConstIter(wordHashSet, allCloudNames, cloudIter)
334  {
335  // Add the name of the cloud(s) to the case file header
336  if (Pstream::master())
337  {
338  ensightCaseFile
339  << (
340  "measured: 1 "
341  + prepend
342  + "****."
343  + cloudIter.key()
344  ).c_str()
345  << nl;
346  }
347 
348  // Create a new hash table for each cloud
349  allCloudFields.insert(cloudIter.key(), HashTable<word>());
350 
351  // Identify the new cloud in the hash table
352  HashTable<HashTable<word>>::iterator newCloudIter =
353  allCloudFields.find(cloudIter.key());
354 
355  // Loop over all times to build list of fields and field types
356  // for each cloud
357  forAll(Times, timeI)
358  {
359  runTime.setTime(Times[timeI], timeI);
360 
361  IOobjectList cloudObjs
362  (
363  mesh,
364  runTime.name(),
365  cloud::prefix/cloudIter.key()
366  );
367 
368  forAllConstIter(IOobjectList, cloudObjs, fieldIter)
369  {
370  const IOobject obj = *fieldIter();
371 
372  if (obj.name() != "positions")
373  {
374  // Add field and field type
375  newCloudIter().insert
376  (
377  obj.name(),
378  obj.headerClassName()
379  );
380  }
381  }
382  }
383  }
384 
385  label nTimeSteps = 0;
386  forAll(Times, timeIndex)
387  {
388  nTimeSteps++;
389  runTime.setTime(Times[timeIndex], timeIndex);
390 
392  word timeFile = prepend + timeName;
393 
394  Info<< "Translating time = " << runTime.name() << nl;
395 
396  polyMesh::readUpdateState meshState = mesh.readUpdate();
397  if (timeIndex != 0 && meshSubsetter.hasSubMesh())
398  {
399  Info<< "Converting cellZone " << cellZoneName
400  << " only (puts outside faces into patch "
401  << mesh.boundaryMesh()[0].name()
402  << ")" << endl;
403  const cellZone& cz = mesh.cellZones()[cellZoneName];
404  cellSet c0(mesh, "c0", labelHashSet(cz));
405  meshSubsetter.setLargeCellSubset(c0, 0);
406  }
407 
408 
409  if (meshState != polyMesh::UNCHANGED)
410  {
411  eMesh.correct();
412  }
413 
414  if (timeIndex == 0 || meshMoving)
415  {
416  eMesh.write
417  (
418  ensightDir,
419  prepend,
420  timeIndex,
421  meshMoving,
422  ensightCaseFile
423  );
424  }
425 
426 
427  // Start of field data output
428  // ~~~~~~~~~~~~~~~~~~~~~~~~~~
429 
430  if (timeIndex == 0 && Pstream::master())
431  {
432  ensightCaseFile<< nl << "VARIABLE" << nl;
433  }
434 
435 
436  // Cell field data output
437  // ~~~~~~~~~~~~~~~~~~~~~~
438 
439  for (label i=0; i<nVolFieldTypes; i++)
440  {
441  wordList fieldNames = objects.names(volFieldTypes[i]);
442 
443  forAll(fieldNames, j)
444  {
445  const word& fieldName = fieldNames[j];
446 
447  // Check if the field has to be exported
448  if (selectedFields)
449  {
450  if (!findStrings(fieldPatterns, fieldName))
451  {
452  continue;
453  }
454  }
455 
456  #include "checkData.H"
457 
458  if (!variableGood)
459  {
460  continue;
461  }
462 
464  (
465  fieldName,
466  mesh.time().name(),
467  mesh,
470  );
471 
472  if (volFieldTypes[i] == volScalarField::typeName)
473  {
474  volScalarField vf(fieldObject, mesh);
475  ensightField<scalar>
476  (
477  volField(meshSubsetter, vf),
478  eMesh,
479  ensightDir,
480  prepend,
481  timeIndex,
482  binary,
483  nodeValues,
484  ensightCaseFile
485  );
486  }
487  else if (volFieldTypes[i] == volVectorField::typeName)
488  {
489  volVectorField vf(fieldObject, mesh);
490  ensightField<vector>
491  (
492  volField(meshSubsetter, vf),
493  eMesh,
494  ensightDir,
495  prepend,
496  timeIndex,
497  binary,
498  nodeValues,
499  ensightCaseFile
500  );
501  }
502  else if (volFieldTypes[i] == volSphericalTensorField::typeName)
503  {
505  ensightField<sphericalTensor>
506  (
507  volField(meshSubsetter, vf),
508  eMesh,
509  ensightDir,
510  prepend,
511  timeIndex,
512  binary,
513  nodeValues,
514  ensightCaseFile
515  );
516  }
517  else if (volFieldTypes[i] == volSymmTensorField::typeName)
518  {
520  ensightField<symmTensor>
521  (
522  volField(meshSubsetter, vf),
523  eMesh,
524  ensightDir,
525  prepend,
526  timeIndex,
527  binary,
528  nodeValues,
529  ensightCaseFile
530  );
531  }
532  else if (volFieldTypes[i] == volTensorField::typeName)
533  {
534  volTensorField vf(fieldObject, mesh);
535  ensightField<tensor>
536  (
537  volField(meshSubsetter, vf),
538  eMesh,
539  ensightDir,
540  prepend,
541  timeIndex,
542  binary,
543  nodeValues,
544  ensightCaseFile
545  );
546  }
547  }
548  }
549 
550 
551  // Cloud field data output
552  // ~~~~~~~~~~~~~~~~~~~~~~~
553 
554  forAllConstIter(HashTable<HashTable<word>>, allCloudFields, cloudIter)
555  {
556  const word& cloudName = cloudIter.key();
557 
558  fileNameList currentCloudDirs = readDir
559  (
560  runTime.timePath()/regionPrefix/cloud::prefix,
562  );
563 
564  bool cloudExists = inFileNameList(currentCloudDirs, cloudName);
566  (
567  mesh,
568  ensightDir,
569  timeFile,
570  cloudName,
571  cloudExists
572  );
573 
574  forAllConstIter(HashTable<word>, cloudIter(), fieldIter)
575  {
576  const word& fieldName = fieldIter.key();
577  const word& fieldType = fieldIter();
578 
580  (
581  fieldName,
582  mesh.time().name(),
584  mesh,
586  );
587 
588  bool fieldExists = fieldObject.headerOk();
589 
590  if (fieldType == scalarIOField::typeName)
591  {
592  ensightCloudField<scalar>
593  (
594  fieldObject,
595  ensightDir,
596  prepend,
597  timeIndex,
598  cloudName,
599  ensightCaseFile,
600  fieldExists
601  );
602  }
603  else if (fieldType == vectorIOField::typeName)
604  {
605  ensightCloudField<vector>
606  (
607  fieldObject,
608  ensightDir,
609  prepend,
610  timeIndex,
611  cloudName,
612  ensightCaseFile,
613  fieldExists
614  );
615  }
616  else
617  {
618  Info<< "Unable to convert field type " << fieldType
619  << " for field " << fieldName << endl;
620  }
621  }
622  }
623  }
624 
625  #include "ensightCaseTail.H"
626 
627  if (Pstream::master())
628  {
629  delete ensightCaseFilePtr;
630  }
631 
632  Info<< "End\n" << endl;
633 
634  return 0;
635 }
636 
637 
638 // ************************************************************************* //
Istream and Ostream manipulators taking arguments.
#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
bool variableGood
Definition: checkData.H:3
bool meshMoving
static const char *const typeName
Definition: Field.H:105
Generic GeometricField class.
A HashTable with keys but without contents.
Definition: HashSet.H:62
bool insert(const Key &key)
Insert a new entry.
Definition: HashSet.H:111
An STL-conforming hash table.
Definition: HashTable.H:127
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
iterator find(const Key &)
Find and return an iterator set at the hashedEntry.
Definition: HashTable.C:142
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
bool headerOk()
Read header of local object without type-checking.
const word & headerClassName() const
Return name of the class name read from header.
Definition: IOobject.H:316
const word & name() const
Return name.
Definition: IOobject.H:310
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Output to file stream.
Definition: OFstream.H:86
T & last()
Return the last element of the list.
Definition: UListI.H:128
static bool master(const label communicator=0)
Am I the master process.
Definition: UPstream.H:423
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
const fileName & globalCaseName() const
Return case name.
Definition: argListI.H:54
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
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
const fileName & rootPath() const
Return root path.
Definition: argListI.H:42
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:120
A collection of cell labels.
Definition: cellSet.H:51
A subset of mesh cells.
Definition: cellZone.H:64
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:63
A class for handling file names.
Definition: fileName.H:82
Post-processing mesh subset tool. Given the original mesh and the list of selected cells,...
Definition: fvMeshSubset.H:74
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:268
readUpdateState
Enumeration defining the state of the mesh after a read update.
Definition: polyMesh.H:91
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
Foam::word regionName
Foam::tmp< Foam::VolField< Type > > volField(const Foam::fvMeshSubset &, const Foam::VolField< Type > &vf)
Wrapper to get hold of the field or the subsetted field.
int main(int argc, char *argv[])
Definition: financialFoam.C:44
IOobject fieldObject(fieldNames[var2field[nVar]], runTime.name(), mesh, IOobject::MUST_READ, IOobject::NO_WRITE)
static List< word > fieldNames
Definition: globalFoam.H:46
word timeName
Definition: getTimeIndex.H:3
Namespace for OpenFOAM.
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
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
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
void ensightParticlePositions(const polyMesh &mesh, const fileName &dataDir, const fileName &subDir, const word &cloudName, IOstream::streamFormat format)
bool rmDir(const fileName &)
Remove a directory and its contents.
Definition: POSIX.C:1047
word itoa(const label n)
bool isDir(const fileName &, const bool followLink=true)
Does the name exist as a directory in the file system?
Definition: POSIX.C:539
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
HashSet< label, Hash< label > > labelHashSet
A HashSet with label keys.
Definition: HashSet.H:211
static const char nl
Definition: Ostream.H:260
fileNameList readDir(const fileName &, const fileType=fileType::file, const bool filterVariants=true, const bool followLink=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:662
label timeIndex
Definition: getTimeIndex.H:4
objects
Foam::argList args(argc, argv)
const word cloudName(propsDict.lookup("cloudName"))