foamFormatConvert.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  foamFormatConvert
26 
27 Description
28  Converts all IOobjects associated with a case into the format specified
29  in the controlDict.
30 
31  Mainly used to convert binary mesh/field files to ASCII.
32 
33  Problem: any zero-size List written binary gets written as '0'. When
34  reading the file as a dictionary this is interpreted as a label. This
35  is (usually) not a problem when doing patch fields since these get the
36  'uniform', 'nonuniform' prefix. However zone contents are labelLists
37  not labelFields and these go wrong. For now hacked a solution where
38  we detect the keywords in zones and redo the dictionary entries
39  to be labelLists.
40 
41 \*---------------------------------------------------------------------------*/
42 
43 #include "argList.H"
44 #include "timeSelector.H"
45 #include "Time.H"
46 #include "volFields.H"
47 #include "surfaceFields.H"
48 #include "pointFields.H"
49 #include "cellIOList.H"
50 #include "IOobjectList.H"
51 #include "IOPtrList.H"
52 #include "cloud.H"
53 #include "labelIOField.H"
54 #include "scalarIOField.H"
55 #include "sphericalTensorIOField.H"
56 #include "symmTensorIOField.H"
57 #include "tensorIOField.H"
58 #include "labelFieldIOField.H"
59 #include "vectorFieldIOField.H"
60 #include "Cloud.H"
61 #include "passiveParticle.H"
62 #include "fieldDictionary.H"
63 
64 #include "writeMeshObject.H"
65 
66 using namespace Foam;
67 
68 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
69 
70 namespace Foam
71 {
73 }
74 
75 
76 // Hack to do zones which have Lists in them. See above.
77 bool writeZones(const word& name, const fileName& meshDir, Time& runTime)
78 {
79  IOobject io
80  (
81  name,
82  runTime.timeName(),
83  meshDir,
84  runTime,
87  false
88  );
89 
90  bool writeOk = false;
91 
92  if (io.headerOk())
93  {
94  Info<< " Reading " << io.headerClassName()
95  << " : " << name << endl;
96 
97  // Switch off type checking (for reading e.g. faceZones as
98  // generic list of dictionaries).
99  const word oldTypeName = IOPtrList<entry>::typeName;
101 
103 
104  forAll(meshObject, i)
105  {
106  if (meshObject[i].isDict())
107  {
108  dictionary& d = meshObject[i].dict();
109 
110  if (d.found("faceLabels"))
111  {
112  d.set("faceLabels", labelList(d.lookup("faceLabels")));
113  }
114 
115  if (d.found("flipMap"))
116  {
117  d.set("flipMap", boolList(d.lookup("flipMap")));
118  }
119 
120  if (d.found("cellLabels"))
121  {
122  d.set("cellLabels", labelList(d.lookup("cellLabels")));
123  }
124 
125  if (d.found("pointLabels"))
126  {
127  d.set("pointLabels", labelList(d.lookup("pointLabels")));
128  }
129  }
130  }
131 
132  const_cast<word&>(IOPtrList<entry>::typeName) = oldTypeName;
133  // Fake type back to what was in field
134  const_cast<word&>(meshObject.type()) = io.headerClassName();
135 
136  Info<< " Writing " << name << endl;
137 
138  // Force writing as ascii
139  writeOk = meshObject.regIOobject::writeObject
140  (
143  runTime.writeCompression(),
144  true
145  );
146  }
147 
148  return writeOk;
149 }
150 
151 
152 // Reduction for non-empty strings
153 class uniqueEqOp
154 {
155  public:
156  void operator()(stringList& x, const stringList& y) const
157  {
158  stringList newX(x.size()+y.size());
159  label n = 0;
160  forAll(x, i)
161  {
162  if (!x[i].empty())
163  {
164  newX[n++] = x[i];
165  }
166  }
167  forAll(y, i)
168  {
169  if (!y[i].empty() && findIndex(x, y[i]) == -1)
170  {
171  newX[n++] = y[i];
172  }
173  }
174  newX.setSize(n);
175  x.transfer(newX);
176  }
177 };
178 
179 
180 template<class T>
181 bool writeOptionalMeshObject
182 (
183  const word& name,
184  const fileName& meshDir,
185  Time& runTime,
186  const bool write
187 )
188 {
189  IOobject io
190  (
191  name,
192  runTime.timeName(),
193  meshDir,
194  runTime,
197  false
198  );
199 
200  bool writeOk = false;
201 
202  bool haveFile = io.headerOk();
203 
204  // Make sure all know if there is a valid class name
205  stringList classNames(1, io.headerClassName());
206  combineReduce(classNames, uniqueEqOp());
207 
208  // Check for correct type
209  if (classNames[0] == T::typeName)
210  {
211  Info<< " Reading " << classNames[0]
212  << " : " << name << endl;
213  T meshObject(io, write && haveFile);
214 
215  Info<< " Writing " << name << endl;
216  writeOk = meshObject.regIOobject::write(write && haveFile);
217  }
218 
219  return writeOk;
220 }
221 
222 
223 int main(int argc, char *argv[])
224 {
227  (
228  "noConstant",
229  "exclude the 'constant/' dir in the times list"
230  );
231 
232  #include "addRegionOption.H"
233  #include "setRootCase.H"
234 
235  // enable noConstant by switching
236  if (!args.optionFound("noConstant"))
237  {
238  args.setOption("constant", "");
239  }
240  else
241  {
242  args.unsetOption("constant");
243  Info<< "Excluding the constant directory." << nl << endl;
244  }
245 
246 
247  #include "createTime.H"
248  // Optional mesh (used to read Clouds)
250 
251 
252  // Make sure we do not use the master-only reading since we read
253  // fields (different per processor) as dictionaries.
255 
256 
257  fileName meshDir = polyMesh::meshSubDir;
258  fileName regionPrefix = "";
259  word regionName = polyMesh::defaultRegion;
260  if (args.optionReadIfPresent("region", regionName))
261  {
262  Info<< "Using region " << regionName << nl << endl;
263  regionPrefix = regionName;
264  meshDir = regionName/polyMesh::meshSubDir;
265  }
266 
268 
269  forAll(timeDirs, timeI)
270  {
271  runTime.setTime(timeDirs[timeI], timeI);
272  Info<< "Time = " << runTime.userTimeName() << endl;
273 
274  // Convert all the standard mesh files
275  writeMeshObject<cellCompactIOList, cellIOList>
276  (
277  "cells",
278  meshDir,
279  runTime
280  );
281  writeMeshObject<labelIOList>("owner", meshDir, runTime);
282  writeMeshObject<labelIOList>("neighbour", meshDir, runTime);
283  writeMeshObject<faceCompactIOList, faceIOList>
284  (
285  "faces",
286  meshDir,
287  runTime
288  );
289  writeMeshObject<pointIOField>("points", meshDir, runTime);
290  // Write boundary in ascii. This is only needed for fileHandler to
291  // kick in. Should not give problems since always writing ascii.
292  writeZones("boundary", meshDir, runTime);
293  writeMeshObject<labelIOList>("pointProcAddressing", meshDir, runTime);
294  writeMeshObject<labelIOList>("faceProcAddressing", meshDir, runTime);
295  writeMeshObject<labelIOList>("cellProcAddressing", meshDir, runTime);
296 
297  // foamyHexMesh vertices
298  writeMeshObject<pointIOField>
299  (
300  "internalDelaunayVertices",
301  regionPrefix,
302  runTime
303  );
304 
305  if (runTime.writeFormat() == IOstream::ASCII)
306  {
307  // Only do zones when converting from binary to ascii
308  // The other way gives problems since working on dictionary level.
309  writeZones("cellZones", meshDir, runTime);
310  writeZones("faceZones", meshDir, runTime);
311  writeZones("pointZones", meshDir, runTime);
312  }
313 
314  // Get list of objects from the database
315  IOobjectList objects(runTime, runTime.timeName(), regionPrefix);
316 
318  {
319  const word& headerClassName = iter()->headerClassName();
320 
321  if
322  (
323  headerClassName == volScalarField::typeName
324  || headerClassName == volVectorField::typeName
325  || headerClassName == volSphericalTensorField::typeName
326  || headerClassName == volSymmTensorField::typeName
327  || headerClassName == volTensorField::typeName
328 
329  || headerClassName == surfaceScalarField::typeName
330  || headerClassName == surfaceVectorField::typeName
331  || headerClassName == surfaceSphericalTensorField::typeName
332  || headerClassName == surfaceSymmTensorField::typeName
333  || headerClassName == surfaceTensorField::typeName
334 
335  || headerClassName == pointScalarField::typeName
336  || headerClassName == pointVectorField::typeName
337  || headerClassName == pointSphericalTensorField::typeName
338  || headerClassName == pointSymmTensorField::typeName
339  || headerClassName == pointTensorField::typeName
340 
341  || headerClassName == volScalarField::Internal::typeName
342  || headerClassName == volVectorField::Internal::typeName
343  || headerClassName == volSphericalTensorField::Internal::typeName
344  || headerClassName == volSymmTensorField::Internal::typeName
345  || headerClassName == volTensorField::Internal::typeName
346  )
347  {
348  Info<< " Reading " << headerClassName
349  << " : " << iter()->name() << endl;
350 
351  fieldDictionary fDict(*iter(), headerClassName);
352 
353  Info<< " Writing " << iter()->name() << endl;
354  fDict.regIOobject::write();
355  }
356  }
357 
358 
359 
360  // Check for lagrangian
361  const fileName lagrangianDir
362  (
363  fileHandler().filePath
364  (
365  runTime.timePath()
366  / regionPrefix
367  / cloud::prefix
368  )
369  );
370  stringList lagrangianDirs
371  (
372  lagrangianDir == fileName::null ? 0 : 1,
373  lagrangianDir
374  );
375 
376  combineReduce(lagrangianDirs, uniqueEqOp());
377 
378  if (!lagrangianDirs.empty())
379  {
380  if (meshPtr.valid())
381  {
382  meshPtr().readUpdate();
383  }
384  else
385  {
386  Info<< " Create polyMesh for time = "
387  << runTime.timeName() << endl;
388 
389  meshPtr.reset
390  (
391  new polyMesh
392  (
393  IOobject
394  (
395  regionName,
396  runTime.timeName(),
397  runTime,
399  )
400  )
401  );
402  }
403 
404  stringList cloudDirs
405  (
407  (
408  lagrangianDirs[0],
410  )
411  );
412 
413  combineReduce(cloudDirs, uniqueEqOp());
414 
415  forAll(cloudDirs, i)
416  {
417  fileName dir(cloud::prefix/cloudDirs[i]);
418 
419  Cloud<passiveParticle> parcels(meshPtr(), cloudDirs[i], false);
420 
421  parcels.writeObject
422  (
423  runTime.writeFormat(),
425  runTime.writeCompression(),
426  parcels.size()
427  );
428 
429 
430  // Do local scan for valid cloud objects
431  IOobjectList sprayObjs(runTime, runTime.timeName(), dir);
432 
433  // Combine with all other cloud objects
434  stringList sprayFields(sprayObjs.sortedToc());
435  combineReduce(sprayFields, uniqueEqOp());
436 
437  forAll(sprayFields, fieldi)
438  {
439  const word& name = sprayFields[fieldi];
440 
441  // Note: try the various field types. Make sure to
442  // exit once successful conversion to avoid re-read
443  // converted file.
444 
445  if
446  (
447  name == "positions"
448  || name == "origProcId"
449  || name == "origId"
450  )
451  {
452  continue;
453  }
454 
455  bool writeOk = writeOptionalMeshObject<labelIOField>
456  (
457  name,
458  dir,
459  runTime,
460  parcels.size() > 0
461  );
462  if (writeOk) continue;
463 
464  writeOk = writeOptionalMeshObject<scalarIOField>
465  (
466  name,
467  dir,
468  runTime,
469  parcels.size() > 0
470  );
471  if (writeOk) continue;
472 
473  writeOk = writeOptionalMeshObject<vectorIOField>
474  (
475  name,
476  dir,
477  runTime,
478  parcels.size() > 0
479  );
480  if (writeOk) continue;
481 
482  writeOk = writeOptionalMeshObject<sphericalTensorIOField>
483  (
484  name,
485  dir,
486  runTime,
487  parcels.size() > 0
488  );
489  if (writeOk) continue;
490 
491  writeOk = writeOptionalMeshObject<symmTensorIOField>
492  (
493  name,
494  dir,
495  runTime,
496  parcels.size() > 0
497  );
498  if (writeOk) continue;
499 
500  writeOk = writeOptionalMeshObject<tensorIOField>
501  (
502  name,
503  dir,
504  runTime,
505  parcels.size() > 0
506  );
507  if (writeOk) continue;
508 
509  writeOk = writeOptionalMeshObject<labelFieldIOField>
510  (
511  name,
512  dir,
513  runTime,
514  parcels.size() > 0
515  );
516  if (writeOk) continue;
517 
518  writeOk = writeOptionalMeshObject<vectorFieldIOField>
519  (
520  name,
521  dir,
522  runTime,
523  parcels.size() > 0
524  );
525 
526  if (!writeOk)
527  {
528  Info<< " Failed converting " << name << endl;
529  }
530  }
531  }
532  }
533 
534  Info<< endl;
535  }
536 
537  Info<< "End\n" << endl;
538 
539  return 0;
540 }
541 
542 
543 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:663
bool setOption(const word &opt, const string &param="")
Set option directly (use with caution)
Definition: argList.C:1112
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void combineReduce(const List< UPstream::commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
Foam::word regionName
A class for handling file names.
Definition: fileName.H:79
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
void reset(T *=nullptr)
If object pointer already set, delete object and set to given.
Definition: autoPtrI.H:114
bool unsetOption(const word &opt)
Unset option directly (use with caution)
Definition: argList.C:1188
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:156
Write a mesh object in the format specified in controlDict.
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
static const char *const typeName
Definition: Field.H:105
static fvMesh * meshPtr
Definition: globalFoam.H:52
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:328
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
static const fileName null
An empty fileName.
Definition: fileName.H:97
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:325
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:114
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:204
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
static word timeName(const scalar, const int precision=curPrecision_)
Return time name of given scalar time.
Definition: Time.C:666
A class for handling words, derived from string.
Definition: word.H:59
static const word null
An empty word.
Definition: word.H:77
List< label > labelList
A List of labels.
Definition: labelList.H:56
const fileOperation & fileHandler()
Get current file handler.
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
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
defineTemplateTypeNameAndDebug(IOPtrList< ensightPart >, 0)
Read field as dictionary (without mesh).
static const char nl
Definition: Ostream.H:260
objects
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:231
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:62
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
readUpdateState readUpdate(const stitchType stitch=stitchType::geometric)
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:713
void setSize(const label)
Reset size of List.
Definition: List.C:281
A PtrList of objects of type <T> with automated input and output.
Definition: IOPtrList.H:50
static const versionNumber currentVersion
Current version number.
Definition: IOstream.H:203
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1291
messageStream Info
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:76
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:118
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
void transfer(List< T > &)
Transfer the contents of the argument List into this list.
Definition: List.C:342
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
Namespace for OpenFOAM.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:864