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-2021 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  (
163  IOstreamOption(IOstream::ASCII, compression),
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  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 
300  fileName meshDir = polyMesh::meshSubDir;
301  fileName regionPrefix = "";
303  if (args.readIfPresent("region", regionName))
304  {
305  Info<< "Using region " << regionName << nl << endl;
306  regionPrefix = regionName;
308  }
309 
311 
312  forAll(timeDirs, timeI)
313  {
314  runTime.setTime(timeDirs[timeI], timeI);
315  Info<< "Time = " << runTime.timeName() << endl;
316 
317  // Convert all the standard mesh files
318  writeMeshObject<cellCompactIOList, cellIOList>
319  (
320  "cells",
321  meshDir,
322  runTime,
323  false // do not check typeName since varies between binary/ascii
324  );
325  writeMeshObject<labelIOList>("owner", meshDir, runTime);
326  writeMeshObject<labelIOList>("neighbour", meshDir, runTime);
327  writeMeshObject<faceCompactIOList, faceIOList>
328  (
329  "faces",
330  meshDir,
331  runTime,
332  false // do not check typeName since varies between binary/ascii
333  );
334  writeMeshObject<pointIOField>("points", meshDir, runTime);
335  // Write boundary in ascii. This is only needed for fileHandler to
336  // kick in. Should not give problems since always writing ascii.
337  writeZones("boundary", meshDir, runTime, IOstreamOption::UNCOMPRESSED);
338  writeMeshObject<labelIOList>("pointProcAddressing", meshDir, runTime);
339  writeMeshObject<labelIOList>("faceProcAddressing", meshDir, runTime);
340  writeMeshObject<labelIOList>("cellProcAddressing", meshDir, runTime);
341  writeMeshObject<labelIOList>
342  (
343  "boundaryProcAddressing",
344  meshDir,
345  runTime
346  );
347 
348  // foamyHexMesh vertices
349  writeMeshObject<pointIOField>
350  (
351  "internalDelaunayVertices",
352  regionPrefix,
353  runTime
354  );
355 
357  {
358  // Only do zones when converting from binary to ascii
359  // The other way gives problems since working on dictionary level.
360  const IOstreamOption::compressionType compress =
362  writeZones("cellZones", meshDir, runTime, compress);
363  writeZones("faceZones", meshDir, runTime, compress);
364  writeZones("pointZones", meshDir, runTime, compress);
365  }
366 
367  // Get list of objects from the database
368  IOobjectList objects(runTime, runTime.timeName(), regionPrefix);
369 
370  forAllConstIters(objects, iter)
371  {
372  const word& headerClassName = (*iter)->headerClassName();
373 
374  if
375  (
376  headerClassName == volScalarField::typeName
377  || headerClassName == volVectorField::typeName
378  || headerClassName == volSphericalTensorField::typeName
379  || headerClassName == volSymmTensorField::typeName
380  || headerClassName == volTensorField::typeName
381 
382  || headerClassName == surfaceScalarField::typeName
383  || headerClassName == surfaceVectorField::typeName
384  || headerClassName == surfaceSphericalTensorField::typeName
385  || headerClassName == surfaceSymmTensorField::typeName
386  || headerClassName == surfaceTensorField::typeName
387 
388  || headerClassName == pointScalarField::typeName
389  || headerClassName == pointVectorField::typeName
390  || headerClassName == pointSphericalTensorField::typeName
391  || headerClassName == pointSymmTensorField::typeName
392  || headerClassName == pointTensorField::typeName
393 
394  || headerClassName == volScalarField::Internal::typeName
395  || headerClassName == volVectorField::Internal::typeName
396  || headerClassName == volSphericalTensorField::Internal::typeName
397  || headerClassName == volSymmTensorField::Internal::typeName
398  || headerClassName == volTensorField::Internal::typeName
399  )
400  {
401  Info<< " Reading " << headerClassName
402  << " : " << (*iter)->name() << endl;
403 
404  fieldDictionary fDict(*iter(), headerClassName);
405 
406  Info<< " Writing " << (*iter)->name() << endl;
407  fDict.regIOobject::write();
408  }
409  }
410 
411 
412 
413  // Check for lagrangian
414  fileNameList lagrangianDirs
415  (
416  1,
417  fileHandler().filePath
418  (
419  runTime.timePath()
420  / regionPrefix
421  / cloud::prefix
422  )
423  );
424 
425  combineReduce(lagrangianDirs, uniqueEqOp<fileName>());
426 
427  if (!lagrangianDirs.empty())
428  {
429  if (meshPtr)
430  {
431  meshPtr->readUpdate();
432  }
433  else
434  {
435  Info<< " Create polyMesh for time = "
436  << runTime.timeName() << endl;
437 
438  meshPtr.reset
439  (
440  new polyMesh
441  (
442  IOobject
443  (
445  runTime.timeName(),
446  runTime,
448  )
449  )
450  );
451  }
452 
453  fileNameList cloudDirs
454  (
456  (
457  lagrangianDirs[0],
459  )
460  );
461 
462  combineReduce(cloudDirs, uniqueEqOp<fileName>());
463 
464  forAll(cloudDirs, i)
465  {
466  fileName dir(cloud::prefix/cloudDirs[i]);
467 
468  Cloud<passiveParticle> parcels(meshPtr(), cloudDirs[i], false);
469 
470  parcels.writeObject
471  (
473  (
476  ),
477  parcels.size()
478  );
479 
480 
481  // Do local scan for valid cloud objects
482  wordList cloudFields
483  (
485  );
486 
487  // Combine with all other cloud objects
488  combineReduce(cloudFields, uniqueEqOp<word>());
489 
490  for (const word& name : cloudFields)
491  {
492  // Note: try the various field types. Make sure to
493  // exit once successful conversion to avoid re-read
494  // converted file.
495 
496  if
497  (
498  name == "positions"
499  || name == "coordinates"
500  || name == "origProcId"
501  || name == "origId"
502  )
503  {
504  continue;
505  }
506 
507  bool writeOk = writeOptionalMeshObject<labelIOField>
508  (
509  name,
510  dir,
511  runTime,
512  parcels.size() > 0
513  );
514  if (writeOk) continue;
515 
516  writeOk = writeOptionalMeshObject<scalarIOField>
517  (
518  name,
519  dir,
520  runTime,
521  parcels.size() > 0
522  );
523  if (writeOk) continue;
524 
525  writeOk = writeOptionalMeshObject<vectorIOField>
526  (
527  name,
528  dir,
529  runTime,
530  parcels.size() > 0
531  );
532  if (writeOk) continue;
533 
534  writeOk = writeOptionalMeshObject<sphericalTensorIOField>
535  (
536  name,
537  dir,
538  runTime,
539  parcels.size() > 0
540  );
541  if (writeOk) continue;
542 
543  writeOk = writeOptionalMeshObject<symmTensorIOField>
544  (
545  name,
546  dir,
547  runTime,
548  parcels.size() > 0
549  );
550  if (writeOk) continue;
551 
552  writeOk = writeOptionalMeshObject<tensorIOField>
553  (
554  name,
555  dir,
556  runTime,
557  parcels.size() > 0
558  );
559  if (writeOk) continue;
560 
561  writeOk = writeOptionalMeshObject<labelFieldIOField>
562  (
563  name,
564  dir,
565  runTime,
566  parcels.size() > 0
567  );
568  if (writeOk) continue;
569 
570  writeOk = writeOptionalMeshObject<vectorFieldIOField>
571  (
572  name,
573  dir,
574  runTime,
575  parcels.size() > 0
576  );
577 
578  if (!writeOk)
579  {
580  Info<< " Failed converting " << name << endl;
581  }
582  }
583  }
584  }
585 
586  Info<< endl;
587  }
588 
590 
591  Info<< "End\n" << endl;
592 
593  return 0;
594 }
595 
596 
597 // ************************************************************************* //
Foam::IOstreamOption::UNCOMPRESSED
@ UNCOMPRESSED
compression = false
Definition: IOstreamOption.H:75
Foam::Time::writeCompression
IOstream::compressionType writeCompression() const
Definition: Time.H:389
Foam::IOobject::NO_WRITE
@ NO_WRITE
Definition: IOobject.H:191
Foam::labelList
List< label > labelList
A List of labels.
Definition: List.H:63
volFields.H
Foam::HashTable::size
label size() const noexcept
Definition: HashTableI.H:45
runTime
engineTime & runTime
Definition: createEngineTime.H:13
Foam::cloud::prefix
static const word prefix
Definition: cloud.H:83
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:165
Foam::Time
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
Foam::IOobject::timeStamp
@ timeStamp
Definition: IOobject.H:197
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:63
Foam::fileName
A class for handling file names.
Definition: fileName.H:71
Foam::IOField
A primitive field of type <T> with automated input and output.
Definition: foamVtkLagrangianWriter.H:57
Foam::IOPtrList
A PtrList of objects of type <T> with automated input and output.
Definition: IOPtrList.H:49
Foam::polyMesh::defaultRegion
static word defaultRegion
Definition: polyMesh.H:314
Foam::fieldDictionary
Read field as dictionary (without mesh).
Definition: fieldDictionary.H:43
Foam::dictionary::found
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionaryI.H:80
Foam::Time::writeFormat
IOstream::streamFormat writeFormat() const
Definition: Time.H:383
Foam::polyMesh::meshSubDir
static word meshSubDir
Definition: polyMesh.H:317
Foam::Time::timeName
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:773
Foam::IOobject::fileModificationChecking
static fileCheckTypes fileModificationChecking
Definition: IOobject.H:299
Cloud.H
cloud.H
tensorIOField.H
Foam::argList::addNote
static void addNote(const string &note)
Definition: argList.C:405
Foam::IOobjectList::sortedNames
wordList sortedNames() const
Definition: IOobjectList.C:338
Foam::combineReduce
void combineReduce(const List< UPstream::commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
Definition: PstreamCombineReduceOps.H:48
Foam::boolList
List< bool > boolList
A List of bools.
Definition: List.H:61
Foam::fileHandler
const fileOperation & fileHandler()
Definition: fileOperation.C:1478
IOobjectList.H
Foam::defineTemplateTypeNameAndDebug
defineTemplateTypeNameAndDebug(faScalarMatrix, 0)
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
surfaceFields.H
Foam::surfaceFields.
scalarIOField.H
Foam::dictionary::set
entry * set(entry *entryPtr)
Definition: dictionary.C:773
writeMeshObject.H
Write a mesh object in the format specified in controlDict.
Foam::argList::readIfPresent
bool readIfPresent(const word &optName, T &val) const
Definition: argListI.H:316
Foam::polyMesh
Mesh consisting of general polyhedral cells.
Definition: polyMesh.H:73
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
n
label n
Definition: TABSMDCalcMethod2.H:31
regionName
Foam::word regionName
Definition: createNamedDynamicFvMesh.H:1
IOPtrList.H
vectorFieldIOField.H
Foam::argList::setOption
bool setOption(const word &optName, const string &param="")
Definition: argList.C:1779
Foam::Info
messageStream Info
Foam::Time::timePath
fileName timePath() const
Definition: Time.H:371
passiveParticle.H
argList.H
Foam::T
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Definition: FieldFieldFunctions.C:51
Foam::IOstreamOption
The IOstreamOption is a simple container for options an IOstream can normally have.
Definition: IOstreamOption.H:59
addRegionOption.H
Foam::ZoneMesh< cellZone, polyMesh >
labelFieldIOField.H
Foam::dictionary::lookup
ITstream & lookup(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Definition: dictionary.C:379
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:119
Foam
Definition: atmBoundaryLayer.C:26
Foam::entry::disableFunctionEntries
static int disableFunctionEntries
Definition: entry.H:123
Foam::IOobjectList
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:51
Foam::argList::unsetOption
bool unsetOption(const word &optName)
Definition: argList.C:1805
Foam::meshObject
Definition: MeshObject.H:134
Foam::autoPtr::reset
void reset(autoPtr< T > &&other) noexcept
Foam::argList::addBoolOption
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Definition: argList.C:317
Time.H
Foam::autoPtr
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Definition: HashPtrTable.H:49
setRootCase.H
Foam::IOstreamOption::ASCII
@ ASCII
"ascii" (normal default)
Definition: IOstreamOption.H:68
Foam::nl
constexpr char nl
Definition: Ostream.H:424
forAllConstIters
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
Foam::timeSelector::addOptions
static void addOptions(const bool constant=true, const bool withZero=false)
Definition: timeSelector.C:95
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: BitOps.H:58
Foam::Time::setTime
virtual void setTime(const Time &t)
Definition: Time.C:996
Foam::Cloud< passiveParticle >
Foam::fileName::DIRECTORY
@ DIRECTORY
A directory.
Definition: fileName.H:82
Foam::word::null
static const word null
Definition: word.H:78
timeSelector.H
createTime.H
x
x
Definition: LISASMDCalcMethod2.H:52
Foam::name
word name(const expressions::valueTypeCode typeCode)
Definition: exprTraits.C:52
symmTensorIOField.H
Foam::IOstreamOption::compressionType
compressionType
Definition: IOstreamOption.H:73
sphericalTensorIOField.H
Foam::DimensionedField::typeName
const word typeName("volScalarField::Internal")
Foam::timeSelector::select0
static instantList select0(Time &runTime, const argList &args)
Definition: timeSelector.C:228
args
Foam::argList args(argc, argv)
Foam::readDir
fileNameList readDir(const fileName &directory, const fileName::Type type=fileName::FILE, const bool filtergz=true, const bool followLink=true)
Definition: POSIX.C:878
meshPtr
Foam::autoPtr< Foam::fvMesh > meshPtr(nullptr)
cellIOList.H
pointFields.H
fieldDictionary.H
labelIOField.H
Foam::argList::found
bool found(const word &optName) const
Definition: argListI.H:171
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:181