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-2024 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.name(),
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 
102  IOPtrList<entry> meshObject(io);
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 
134  // Fake type back to what was in field
135  const_cast<word&>(meshObject.type()) = io.headerClassName();
136 
137  Info<< " Writing " << name << endl;
138 
139  // Force writing as ascii
140  writeOk = meshObject.regIOobject::writeObject
141  (
144  runTime.writeCompression(),
145  true
146  );
147  }
148 
149  return writeOk;
150 }
151 
152 
153 // Reduction for non-empty strings
154 class uniqueEqOp
155 {
156  public:
157  void operator()(stringList& x, const stringList& y) const
158  {
159  stringList newX(x.size()+y.size());
160  label n = 0;
161  forAll(x, i)
162  {
163  if (!x[i].empty())
164  {
165  newX[n++] = x[i];
166  }
167  }
168  forAll(y, i)
169  {
170  if (!y[i].empty() && findIndex(x, y[i]) == -1)
171  {
172  newX[n++] = y[i];
173  }
174  }
175  newX.setSize(n);
176  x.transfer(newX);
177  }
178 };
179 
180 
181 template<class T>
182 bool writeOptionalMeshObject
183 (
184  const word& name,
185  const fileName& meshDir,
186  Time& runTime,
187  const bool write
188 )
189 {
190  IOobject io
191  (
192  name,
193  runTime.name(),
194  meshDir,
195  runTime,
198  false
199  );
200 
201  bool writeOk = false;
202 
203  const bool haveFile = io.headerOk();
204 
205  // Make sure all know if there is a valid class name
206  stringList classNames(1, io.headerClassName());
207  combineReduce(classNames, uniqueEqOp());
208 
209  // Check for correct type
210  if (classNames[0] == T::typeName)
211  {
212  Info<< " Reading " << classNames[0]
213  << " : " << name << endl;
214  T meshObject(io, write && haveFile);
215 
216  Info<< " Writing " << name << endl;
217  writeOk = meshObject.regIOobject::write(write && haveFile);
218  }
219 
220  return writeOk;
221 }
222 
223 
224 int main(int argc, char *argv[])
225 {
228  (
229  "noConstant",
230  "exclude the 'constant/' dir in the times list"
231  );
232 
233  #include "addRegionOption.H"
234  #include "setRootCase.H"
235 
236  // enable noConstant by switching
237  if (!args.optionFound("noConstant"))
238  {
239  args.setOption("constant", "");
240  }
241  else
242  {
243  args.unsetOption("constant");
244  Info<< "Excluding the constant directory." << nl << endl;
245  }
246 
247  #include "createTime.H"
248 
249  // Optional mesh (used to read Clouds)
251 
252  // Make sure we do not use the master-only reading since we read
253  // fields (different per processor) as dictionaries.
255 
256  fileName meshDir = polyMesh::meshSubDir;
257  fileName regionPrefix = "";
259  if (args.optionReadIfPresent("region", regionName))
260  {
261  Info<< "Using region " << regionName << nl << endl;
262  regionPrefix = regionName;
264  }
265 
267 
268  forAll(timeDirs, timeI)
269  {
270  runTime.setTime(timeDirs[timeI], timeI);
271  Info<< "Time = " << runTime.userTimeName() << endl;
272 
273  // Convert all the standard mesh files
274  writeMeshObject<cellCompactIOList, cellIOList>
275  (
276  "cells",
277  meshDir,
278  runTime
279  );
280  writeMeshObject<labelIOList>("owner", meshDir, runTime);
281  writeMeshObject<labelIOList>("neighbour", meshDir, runTime);
282  if (!writeMeshObject<faceCompactIOList>("faces", meshDir, runTime))
283  {
284  writeMeshObject<faceCompactIOList, faceIOList>
285  (
286  "faces",
287  meshDir,
288  runTime
289  );
290  }
291  writeMeshObject<pointIOField>("points", meshDir, runTime);
292  // Write boundary in ascii. This is only needed for fileHandler to
293  // kick in. Should not give problems since always writing ascii.
294  writeZones("boundary", meshDir, runTime);
295  writeMeshObject<labelIOList>("pointProcAddressing", meshDir, runTime);
296  writeMeshObject<labelIOList>("faceProcAddressing", meshDir, runTime);
297  writeMeshObject<labelIOList>("cellProcAddressing", meshDir, runTime);
298 
299  if (runTime.writeFormat() == IOstream::ASCII)
300  {
301  // Only do zones when converting from binary to ascii
302  // The other way gives problems since working on dictionary level.
303  writeZones("cellZones", meshDir, runTime);
304  writeZones("faceZones", meshDir, runTime);
305  writeZones("pointZones", meshDir, runTime);
306  }
307 
308  // Get list of objects from the database
309  IOobjectList objects(runTime, runTime.name(), regionPrefix);
310 
312  {
313  const word& headerClassName = iter()->headerClassName();
314 
315  if
316  (
317  headerClassName == volScalarField::typeName
318  || headerClassName == volVectorField::typeName
319  || headerClassName == volSphericalTensorField::typeName
320  || headerClassName == volSymmTensorField::typeName
321  || headerClassName == volTensorField::typeName
322 
323  || headerClassName == surfaceScalarField::typeName
324  || headerClassName == surfaceVectorField::typeName
325  || headerClassName == surfaceSphericalTensorField::typeName
326  || headerClassName == surfaceSymmTensorField::typeName
327  || headerClassName == surfaceTensorField::typeName
328 
329  || headerClassName == pointScalarField::typeName
330  || headerClassName == pointVectorField::typeName
331  || headerClassName == pointSphericalTensorField::typeName
332  || headerClassName == pointSymmTensorField::typeName
333  || headerClassName == pointTensorField::typeName
334 
335  || headerClassName == volScalarField::Internal::typeName
336  || headerClassName == volVectorField::Internal::typeName
337  || headerClassName == volSphericalTensorField::Internal::typeName
338  || headerClassName == volSymmTensorField::Internal::typeName
339  || headerClassName == volTensorField::Internal::typeName
340  )
341  {
342  Info<< " Reading " << headerClassName
343  << " : " << iter()->name() << endl;
344 
345  fieldDictionary fDict(*iter(), headerClassName);
346 
347  Info<< " Writing " << iter()->name() << endl;
348  fDict.regIOobject::write();
349  }
350  }
351 
352  // Check for lagrangian
353  const fileName lagrangianDir
354  (
355  fileHandler().filePath
356  (
357  runTime.timePath()/regionPrefix/cloud::prefix
358  )
359  );
360  stringList lagrangianDirs
361  (
362  lagrangianDir == fileName::null ? 0 : 1,
363  lagrangianDir
364  );
365 
366  combineReduce(lagrangianDirs, uniqueEqOp());
367 
368  if (!lagrangianDirs.empty())
369  {
370  if (meshPtr.valid())
371  {
372  meshPtr().readUpdate();
373  }
374  else
375  {
376  Info<< " Create polyMesh for time = "
377  << runTime.name() << endl;
378 
379  meshPtr.reset
380  (
381  new polyMesh
382  (
383  IOobject
384  (
385  regionName,
386  runTime.name(),
387  runTime,
389  )
390  )
391  );
392  }
393 
394  stringList cloudDirs
395  (
397  (
398  lagrangianDirs[0],
400  )
401  );
402 
403  combineReduce(cloudDirs, uniqueEqOp());
404 
405  forAll(cloudDirs, i)
406  {
407  const fileName dir(cloud::prefix/cloudDirs[i]);
408 
409  Cloud<passiveParticle> parcels(meshPtr(), cloudDirs[i], false);
410 
411  parcels.writeObject
412  (
413  runTime.writeFormat(),
415  runTime.writeCompression(),
416  parcels.size()
417  );
418 
419 
420  // Do local scan for valid cloud objects
421  IOobjectList sprayObjs(runTime, runTime.name(), dir);
422 
423  // Combine with all other cloud objects
424  stringList sprayFields(sprayObjs.sortedToc());
425  combineReduce(sprayFields, uniqueEqOp());
426 
427  forAll(sprayFields, fieldi)
428  {
429  const word& name = sprayFields[fieldi];
430 
431  // Note: try the various field types. Make sure to
432  // exit once successful conversion to avoid re-read
433  // converted file.
434 
435  if
436  (
437  name == "positions"
438  || name == "origProcId"
439  || name == "origId"
440  )
441  {
442  continue;
443  }
444 
445  bool writeOk = writeOptionalMeshObject<labelIOField>
446  (
447  name,
448  dir,
449  runTime,
450  parcels.size() > 0
451  );
452  if (writeOk) continue;
453 
454  writeOk = writeOptionalMeshObject<scalarIOField>
455  (
456  name,
457  dir,
458  runTime,
459  parcels.size() > 0
460  );
461  if (writeOk) continue;
462 
463  writeOk = writeOptionalMeshObject<vectorIOField>
464  (
465  name,
466  dir,
467  runTime,
468  parcels.size() > 0
469  );
470  if (writeOk) continue;
471 
472  writeOk = writeOptionalMeshObject<sphericalTensorIOField>
473  (
474  name,
475  dir,
476  runTime,
477  parcels.size() > 0
478  );
479  if (writeOk) continue;
480 
481  writeOk = writeOptionalMeshObject<symmTensorIOField>
482  (
483  name,
484  dir,
485  runTime,
486  parcels.size() > 0
487  );
488  if (writeOk) continue;
489 
490  writeOk = writeOptionalMeshObject<tensorIOField>
491  (
492  name,
493  dir,
494  runTime,
495  parcels.size() > 0
496  );
497  if (writeOk) continue;
498 
499  writeOk = writeOptionalMeshObject<labelFieldIOField>
500  (
501  name,
502  dir,
503  runTime,
504  parcels.size() > 0
505  );
506  if (writeOk) continue;
507 
508  writeOk = writeOptionalMeshObject<vectorFieldIOField>
509  (
510  name,
511  dir,
512  runTime,
513  parcels.size() > 0
514  );
515 
516  if (!writeOk)
517  {
518  Info<< " Failed converting " << name << endl;
519  }
520  }
521  }
522  }
523 
524  Info<< endl;
525  }
526 
527  Info<< "End\n" << endl;
528 
529  return 0;
530 }
531 
532 
533 // ************************************************************************* //
scalar y
label n
#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
static const char *const typeName
Definition: Field.H:106
A PtrList of objects of type <Type> with automated input and output.
Definition: IOPtrList.H:53
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
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:226
static const versionNumber currentVersion
Current version number.
Definition: IOstream.H:203
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:76
virtual void setTime(const Time &)
Reset the time and time-index to those of the given time.
Definition: Time.C:985
word userTimeName() const
Return current user time name with units.
Definition: Time.C:856
IOstream::streamFormat writeFormat() const
Default write format.
Definition: Time.H:287
IOstream::compressionType writeCompression() const
Default write compression.
Definition: Time.H:299
fileName timePath() const
Return current time path.
Definition: Time.H:281
bool setOption(const word &opt, const string &param="")
Set option directly (use with caution)
Definition: argList.C:1112
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
bool unsetOption(const word &opt)
Unset option directly (use with caution)
Definition: argList.C:1188
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:63
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:710
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:1152
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:509
const word & name() const
Return const reference to name.
Read field as dictionary (without mesh).
A class for handling file names.
Definition: fileName.H:82
static const fileName null
An empty fileName.
Definition: fileName.H:97
readUpdateState readUpdate(const stitchType stitch=stitchType::geometric)
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:761
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:80
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:267
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:270
void reset(const label nPoints, const label nInternalFaces, const label nFaces, const label nCells)
Reset this primitiveMesh given the primitive array sizes.
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
static const word null
An empty word.
Definition: word.H:77
int main(int argc, char *argv[])
Definition: financialFoam.C:44
static instantList timeDirs
Definition: globalFoam.H:44
static fvMesh * meshPtr
Definition: globalFoam.H:52
void write(std::ostream &os, const bool binary, List< floatScalar > &fField)
Write floats ascii or binary.
Namespace for OpenFOAM.
const fileOperation & fileHandler()
Get current file handler.
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const word & regionName(const solver &region)
Definition: solver.H:209
messageStream Info
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
defineTemplateTypeNameAndDebug(prghPressure, 0)
label findIndex(const ListType &, typename ListType::const_reference, const label start=0)
Find first occurrence of given element and return index,.
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
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
void combineReduce(const List< UPstream::commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
objects
Foam::argList args(argc, argv)
Foam::surfaceFields.
Write a mesh object in the format specified in controlDict.