foamFormatConvert.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-2017 OpenFOAM Foundation
9  Copyright (C) 2016-2022 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  foamFormatConvert
29 
30 Group
31  grpMiscUtilities
32 
33 Description
34  Converts all IOobjects associated with a case into the format specified
35  in the controlDict.
36 
37  Mainly used to convert binary mesh/field files to ASCII.
38 
39  Problem: any zero-size List written binary gets written as '0'. When
40  reading the file as a dictionary this is interpreted as a label. This
41  is (usually) not a problem when doing patch fields since these get the
42  'uniform', 'nonuniform' prefix. However zone contents are labelLists
43  not labelFields and these go wrong. For now hacked a solution where
44  we detect the keywords in zones and redo the dictionary entries
45  to be labelLists.
46 
47 Usage
48 
49  - foamFormatConvert [OPTION]
50 
51  \param -noConstant \n
52  Exclude the constant/ directory from the times list
53 
54  \param -enableFunctionEntries \n
55  By default all dictionary preprocessing of fields is disabled
56 
57 \*---------------------------------------------------------------------------*/
58 
59 #include "argList.H"
60 #include "timeSelector.H"
61 #include "Time.H"
62 #include "volFields.H"
63 #include "surfaceFields.H"
64 #include "pointFields.H"
65 #include "cellIOList.H"
66 #include "IOobjectList.H"
67 #include "IOPtrList.H"
68 #include "cloud.H"
69 #include "labelIOField.H"
70 #include "scalarIOField.H"
71 #include "sphericalTensorIOField.H"
72 #include "symmTensorIOField.H"
73 #include "tensorIOField.H"
74 #include "labelFieldIOField.H"
75 #include "vectorFieldIOField.H"
76 #include "Cloud.H"
77 #include "passiveParticle.H"
78 #include "fieldDictionary.H"
79 
80 #include "writeMeshObject.H"
81 
82 using namespace Foam;
83 
84 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
85 
86 namespace Foam
87 {
89 }
90 
91 
92 // Hack to do zones which have Lists in them. See above.
93 bool writeZones
94 (
95  const word& name,
96  const fileName& meshDir,
97  Time& runTime,
98  const IOstreamOption::compressionType compression
99 )
100 {
101  IOobject io
102  (
103  name,
104  runTime.timeName(),
105  meshDir,
106  runTime,
109  false
110  );
111 
112  bool writeOk = false;
113 
114  if (io.typeHeaderOk<cellZoneMesh>(false))
115  {
116  Info<< " Reading " << io.headerClassName()
117  << " : " << name << endl;
118 
119  // Switch off type checking (for reading e.g. faceZones as
120  // generic list of dictionaries).
121  const word oldTypeName = IOPtrList<entry>::typeName;
123 
125 
126  forAll(meshObject, i)
127  {
128  if (meshObject[i].isDict())
129  {
130  dictionary& d = meshObject[i].dict();
131 
132  if (d.found("faceLabels"))
133  {
134  d.set("faceLabels", labelList(d.lookup("faceLabels")));
135  }
136 
137  if (d.found("flipMap"))
138  {
139  d.set("flipMap", boolList(d.lookup("flipMap")));
140  }
141 
142  if (d.found("cellLabels"))
143  {
144  d.set("cellLabels", labelList(d.lookup("cellLabels")));
145  }
146 
147  if (d.found("pointLabels"))
148  {
149  d.set("pointLabels", labelList(d.lookup("pointLabels")));
150  }
151  }
152  }
153 
154  const_cast<word&>(IOPtrList<entry>::typeName) = oldTypeName;
155  // Fake type back to what was in field
156  const_cast<word&>(meshObject.type()) = io.headerClassName();
157 
158  Info<< " Writing " << name << endl;
159 
160  // Force writing as ASCII
161  writeOk = meshObject.regIOobject::writeObject
162  (
164  true
165  );
166  }
167 
168  return writeOk;
169 }
170 
171 
172 // Reduction for non-empty strings.
173 template<class StringType>
174 struct uniqueEqOp
175 {
176  void operator()(List<StringType>& x, const List<StringType>& y) const
177  {
178  List<StringType> newX(x.size()+y.size());
179  label n = 0;
180  forAll(x, i)
181  {
182  if (!x[i].empty())
183  {
184  newX[n++] = x[i];
185  }
186  }
187  forAll(y, i)
188  {
189  if (!y[i].empty() && !x.found(y[i]))
190  {
191  newX[n++] = y[i];
192  }
193  }
194  newX.setSize(n);
195  x.transfer(newX);
196  }
197 };
198 
199 
200 template<class T>
201 bool writeOptionalMeshObject
202 (
203  const word& name,
204  const fileName& meshDir,
205  Time& runTime,
206  const bool valid
207 )
208 {
209  IOobject io
210  (
211  name,
212  runTime.timeName(),
213  meshDir,
214  runTime,
217  false
218  );
219 
220  bool writeOk = false;
221 
222  bool haveFile = io.typeHeaderOk<IOField<label>>(false);
223 
224  // Make sure all know if there is a valid class name
225  wordList classNames(1, io.headerClassName());
226  Pstream::combineReduce(classNames, uniqueEqOp<word>());
227 
228  // Check for correct type
229  if (classNames[0] == T::typeName)
230  {
231  Info<< " Reading " << classNames[0]
232  << " : " << name << endl;
233  T meshObject(io, valid && haveFile);
234 
235  Info<< " Writing " << name << endl;
236  writeOk = meshObject.regIOobject::write(valid && haveFile);
237  }
238 
239  return writeOk;
240 }
241 
242 
243 int main(int argc, char *argv[])
244 {
246  (
247  "Converts all IOobjects associated with a case into the format"
248  " specified in the controlDict"
249  );
252  (
253  "noConstant",
254  "Exclude the 'constant/' dir in the times list"
255  );
257  (
258  "enableFunctionEntries",
259  "Enable expansion of dictionary directives - #include, #codeStream etc"
260  );
261 
262  #include "addRegionOption.H"
263  #include "setRootCase.H"
264 
265  // enable noConstant by switching
266  if (!args.found("noConstant"))
267  {
268  args.setOption("constant", "");
269  }
270  else
271  {
272  args.unsetOption("constant");
273  Info<< "Excluding the constant directory." << nl << endl;
274  }
275 
276 
277  #include "createTime.H"
278  // Optional mesh (used to read Clouds)
280 
281  const bool enableEntries = args.found("enableFunctionEntries");
282  if (enableEntries)
283  {
284  Info<< "Allowing dictionary preprocessing ('#include', '#codeStream')."
285  << endl;
286  }
287 
288  const int oldFlag = entry::disableFunctionEntries;
289  if (!enableEntries)
290  {
291  // By default disable dictionary expansion for fields
293  }
294 
295  // Make sure we do not use the master-only reading since we read
296  // fields (different per processor) as dictionaries.
298 
299 
301  if (args.readIfPresent("region", regionName))
302  {
303  Info<< "Using region " << regionName << nl << endl;
304  }
305 
306  const fileName meshDir
307  (
309  );
310 
311 
313 
314  forAll(timeDirs, timeI)
315  {
316  runTime.setTime(timeDirs[timeI], timeI);
317  Info<< "Time = " << runTime.timeName() << endl;
318 
319  // Convert all the standard mesh files
320  writeMeshObject<cellCompactIOList, cellIOList>
321  (
322  "cells",
323  meshDir,
324  runTime,
325  false // do not check typeName since varies between binary/ascii
326  );
327  writeMeshObject<labelIOList>("owner", meshDir, runTime);
328  writeMeshObject<labelIOList>("neighbour", meshDir, runTime);
329  writeMeshObject<faceCompactIOList, faceIOList>
330  (
331  "faces",
332  meshDir,
333  runTime,
334  false // do not check typeName since varies between binary/ascii
335  );
336  writeMeshObject<pointIOField>("points", meshDir, runTime);
337  // Write boundary in ascii. This is only needed for fileHandler to
338  // kick in. Should not give problems since always writing ascii.
339  writeZones("boundary", meshDir, runTime, IOstreamOption::UNCOMPRESSED);
340  writeMeshObject<labelIOList>("pointProcAddressing", meshDir, runTime);
341  writeMeshObject<labelIOList>("faceProcAddressing", meshDir, runTime);
342  writeMeshObject<labelIOList>("cellProcAddressing", meshDir, runTime);
343  writeMeshObject<labelIOList>
344  (
345  "boundaryProcAddressing",
346  meshDir,
347  runTime
348  );
349 
350  // foamyHexMesh vertices
351  writeMeshObject<pointIOField>
352  (
353  "internalDelaunayVertices",
355  runTime
356  );
357 
359  {
360  // Only do zones when converting from binary to ascii
361  // The other way gives problems since working on dictionary level.
362  const IOstreamOption::compressionType compress =
364  writeZones("cellZones", meshDir, runTime, compress);
365  writeZones("faceZones", meshDir, runTime, compress);
366  writeZones("pointZones", meshDir, runTime, compress);
367  }
368 
369  // Get list of objects from the database
370  IOobjectList objects
371  (
372  runTime,
373  runTime.timeName(),
375  );
376 
377  forAllConstIters(objects, iter)
378  {
379  const word& headerClassName = (*iter)->headerClassName();
380 
381  if
382  (
383  headerClassName == volScalarField::typeName
384  || headerClassName == volVectorField::typeName
385  || headerClassName == volSphericalTensorField::typeName
386  || headerClassName == volSymmTensorField::typeName
387  || headerClassName == volTensorField::typeName
388 
389  || headerClassName == surfaceScalarField::typeName
390  || headerClassName == surfaceVectorField::typeName
391  || headerClassName == surfaceSphericalTensorField::typeName
392  || headerClassName == surfaceSymmTensorField::typeName
393  || headerClassName == surfaceTensorField::typeName
394 
395  || headerClassName == pointScalarField::typeName
396  || headerClassName == pointVectorField::typeName
397  || headerClassName == pointSphericalTensorField::typeName
398  || headerClassName == pointSymmTensorField::typeName
399  || headerClassName == pointTensorField::typeName
400 
401  || headerClassName == volScalarField::Internal::typeName
402  || headerClassName == volVectorField::Internal::typeName
403  || headerClassName == volSphericalTensorField::Internal::typeName
404  || headerClassName == volSymmTensorField::Internal::typeName
405  || headerClassName == volTensorField::Internal::typeName
406  )
407  {
408  Info<< " Reading " << headerClassName
409  << " : " << (*iter)->name() << endl;
410 
411  fieldDictionary fDict(*iter(), headerClassName);
412 
413  Info<< " Writing " << (*iter)->name() << endl;
414  fDict.regIOobject::write();
415  }
416  }
417 
418 
419 
420  // Check for lagrangian
421  fileNameList lagrangianDirs
422  (
423  1,
424  fileHandler().filePath
425  (
426  runTime.timePath()
428  / cloud::prefix
429  )
430  );
431 
432  Pstream::combineReduce(lagrangianDirs, uniqueEqOp<fileName>());
433 
434  if (!lagrangianDirs.empty())
435  {
436  if (meshPtr)
437  {
438  meshPtr->readUpdate();
439  }
440  else
441  {
442  Info<< " Create polyMesh for time = "
443  << runTime.timeName() << endl;
444 
445  meshPtr.reset
446  (
447  new polyMesh
448  (
449  IOobject
450  (
452  runTime.timeName(),
453  runTime,
455  )
456  )
457  );
458  }
459 
460  fileNameList cloudDirs
461  (
463  (
464  lagrangianDirs[0],
466  )
467  );
468 
469  Pstream::combineReduce(cloudDirs, uniqueEqOp<fileName>());
470 
471  forAll(cloudDirs, i)
472  {
473  fileName dir(cloud::prefix/cloudDirs[i]);
474 
475  Cloud<passiveParticle> parcels(meshPtr(), cloudDirs[i], false);
476 
477  parcels.writeObject
478  (
480  (
483  ),
484  parcels.size()
485  );
486 
487 
488  // Do local scan for valid cloud objects
489  wordList cloudFields
490  (
491  IOobjectList(runTime, runTime.timeName(), dir).sortedNames()
492  );
493 
494  // Combine with all other cloud objects
495  Pstream::combineReduce(cloudFields, uniqueEqOp<word>());
496 
497  for (const word& name : cloudFields)
498  {
499  // Note: try the various field types. Make sure to
500  // exit once successful conversion to avoid re-read
501  // converted file.
502 
503  if
504  (
505  name == "positions"
506  || name == "coordinates"
507  || name == "origProcId"
508  || name == "origId"
509  )
510  {
511  continue;
512  }
513 
514  bool writeOk = writeOptionalMeshObject<labelIOField>
515  (
516  name,
517  dir,
518  runTime,
519  parcels.size() > 0
520  );
521  if (writeOk) continue;
522 
523  writeOk = writeOptionalMeshObject<scalarIOField>
524  (
525  name,
526  dir,
527  runTime,
528  parcels.size() > 0
529  );
530  if (writeOk) continue;
531 
532  writeOk = writeOptionalMeshObject<vectorIOField>
533  (
534  name,
535  dir,
536  runTime,
537  parcels.size() > 0
538  );
539  if (writeOk) continue;
540 
541  writeOk = writeOptionalMeshObject<sphericalTensorIOField>
542  (
543  name,
544  dir,
545  runTime,
546  parcels.size() > 0
547  );
548  if (writeOk) continue;
549 
550  writeOk = writeOptionalMeshObject<symmTensorIOField>
551  (
552  name,
553  dir,
554  runTime,
555  parcels.size() > 0
556  );
557  if (writeOk) continue;
558 
559  writeOk = writeOptionalMeshObject<tensorIOField>
560  (
561  name,
562  dir,
563  runTime,
564  parcels.size() > 0
565  );
566  if (writeOk) continue;
567 
568  writeOk = writeOptionalMeshObject<labelFieldIOField>
569  (
570  name,
571  dir,
572  runTime,
573  parcels.size() > 0
574  );
575  if (writeOk) continue;
576 
577  writeOk = writeOptionalMeshObject<vectorFieldIOField>
578  (
579  name,
580  dir,
581  runTime,
582  parcels.size() > 0
583  );
584 
585  if (!writeOk)
586  {
587  Info<< " Failed converting " << name << endl;
588  }
589  }
590  }
591  }
592 
593  Info<< endl;
594  }
595 
597 
598  Info<< "End\n" << endl;
599 
600  return 0;
601 }
602 
603 
604 // ************************************************************************* //
Foam::surfaceFields.
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return an entry data stream. FatalIOError if not found, or not a stream. ...
Definition: dictionary.C:379
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:453
A class for handling file names.
Definition: fileName.H:71
List of IOobjects with searching and retrieving facilities. Implemented as a HashTable, so the various sorted methods should be used if traversing in parallel.
Definition: IOobjectList.H:55
compressionType
Compression treatment (UNCOMPRESSED | COMPRESSED)
bool unsetOption(const word &optName)
Unset option directly (use with caution)
Definition: argList.C:1900
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:120
Write a mesh object in the format specified in controlDict.
fileName timePath() const
Return current time path.
Definition: Time.H:470
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:841
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:56
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
"ascii" (normal default)
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:402
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:365
autoPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler.
label size() const noexcept
The number of elements in table.
Definition: HashTableI.H:45
A simple container for options an IOstream can normally have.
Ignore writing from objectRegistry::writeObject()
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:100
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:413
scalar y
Foam::word regionName(Foam::polyMesh::defaultRegion)
static int disableFunctionEntries
Enable or disable use of function entries and variable expansions.
Definition: entry.H:139
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, false)
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
Foam::autoPtr< Foam::dynamicFvMesh > meshPtr
A class for handling words, derived from Foam::string.
Definition: word.H:63
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:397
static const word null
An empty word.
Definition: word.H:84
virtual void setTime(const Time &t)
Reset the time and time-index to those of the given time.
Definition: Time.C:977
IOstreamOption::compressionType writeCompression() const noexcept
The write stream compression.
Definition: Time.H:494
static void combineReduce(const List< commsStruct > &comms, T &value, const CombineOp &cop, const int tag=UPstream::msgType(), const label comm=UPstream::worldComm)
Reduce inplace (cf. MPI Allreduce) applying cop to inplace combine value from different processors...
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (uses typeFilePath to find file) and check its info.
Read field as dictionary (without mesh).
static word timeName(const scalar t, const int precision=precision_)
Return time name of given scalar time formatted with the given precision.
Definition: Time.C:760
defineTemplateTypeNameAndDebug(faScalarMatrix, 0)
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:313
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:234
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const word typeName("volScalarField::Internal")
IOstreamOption::streamFormat writeFormat() const noexcept
The write stream format.
Definition: Time.H:486
A PtrList of objects of type <T> with automated input and output.
Definition: IOPtrList.H:49
const word & headerClassName() const noexcept
Return name of the class name read from header.
Definition: IOobjectI.H:168
bool setOption(const word &optName, const string &param="")
Set option directly (use with caution)
Definition: argList.C:1874
messageStream Info
Information stream (stdout output on master, null elsewhere)
label n
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:777
bool readIfPresent(const word &optName, T &val) const
Read a value from the named option if present.
Definition: argListI.H:316
List< label > labelList
A List of labels.
Definition: List.H:62
fileNameList readDir(const fileName &directory, const fileName::Type type=fileName::Type::FILE, const bool filtergz=true, const bool followLink=true)
Read a directory and return the entries as a fileName List.
Definition: POSIX.C:916
Foam::argList args(argc, argv)
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:166
List< bool > boolList
A List of bools.
Definition: List.H:60
A primitive field of type <T> with automated input and output.
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
Namespace for OpenFOAM.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:93