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-2023 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 "passiveParticleCloud.H"
77 #include "fieldDictionary.H"
78 
79 #include "writeMeshObject.H"
80 
81 using namespace Foam;
82 
83 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
84 
85 namespace Foam
86 {
88 }
89 
90 
91 // Hack to do zones which have Lists in them. See above.
92 bool writeZones
93 (
94  const word& name,
95  const fileName& meshDir,
96  Time& runTime,
98 )
99 {
100  IOobject io
101  (
102  name,
103  runTime.timeName(),
104  meshDir,
105  runTime,
109  );
110 
111  bool writeOk = false;
112 
113  if (io.typeHeaderOk<cellZoneMesh>(false))
114  {
115  Info<< " Reading " << io.headerClassName()
116  << " : " << name << endl;
117 
118  // Switch off type checking (for reading e.g. faceZones as
119  // generic list of dictionaries).
120  const word oldTypeName = IOPtrList<entry>::typeName;
122 
124 
125  for (entry& e : meshObject)
126  {
127  if (e.isDict())
128  {
129  dictionary& d = e.dict();
130 
131  if (d.found("faceLabels"))
132  {
133  d.set("faceLabels", labelList(d.lookup("faceLabels")));
134  }
135 
136  if (d.found("flipMap"))
137  {
138  d.set("flipMap", boolList(d.lookup("flipMap")));
139  }
140 
141  if (d.found("cellLabels"))
142  {
143  d.set("cellLabels", labelList(d.lookup("cellLabels")));
144  }
145 
146  if (d.found("pointLabels"))
147  {
148  d.set("pointLabels", labelList(d.lookup("pointLabels")));
149  }
150  }
151  }
152 
153  const_cast<word&>(IOPtrList<entry>::typeName) = oldTypeName;
154  // Fake type back to what was in field
155  const_cast<word&>(meshObject.type()) = io.headerClassName();
156 
157  Info<< " Writing " << name << endl;
158 
159  // Force writing as ASCII
160  writeOk = meshObject.regIOobject::writeObject
161  (
163  true
164  );
165  }
166 
167  return writeOk;
168 }
169 
170 
171 // Reduction for non-empty strings.
172 template<class StringType>
173 struct uniqueEqOp
174 {
175  void operator()(List<StringType>& x, const List<StringType>& y) const
176  {
177  List<StringType> newX(x.size()+y.size());
178  label n = 0;
179  forAll(x, i)
180  {
181  if (!x[i].empty())
182  {
183  newX[n++] = x[i];
184  }
185  }
186  forAll(y, i)
187  {
188  if (!y[i].empty() && !x.found(y[i]))
189  {
190  newX[n++] = y[i];
191  }
192  }
193  newX.setSize(n);
194  x.transfer(newX);
195  }
196 };
197 
198 
199 template<class Type>
200 bool writeOptionalMeshObject
201 (
202  const word& name,
203  const fileName& meshDir,
204  Time& runTime,
205  bool writeOnProc
206 )
207 {
208  IOobject io
209  (
210  name,
211  runTime.timeName(),
212  meshDir,
213  runTime,
217  );
218 
219  // Check if available and the correct type.
220  // Done as typeHeaderOk<regIOobject> + isHeaderClass to ensure
221  // local-only reading and circumvent is_globalIOobject<Type> check
222  // in typeHeaderOk<Type>
223 
224  bool readOnProc =
225  (
226  io.typeHeaderOk<regIOobject>(false)
227  && io.isHeaderClass<Type>()
228  );
229 
230  bool writeOk = false;
231 
232  if (returnReduceOr(readOnProc))
233  {
234  Info<< " Reading " << Type::typeName << " : " << name << endl;
235  Type meshObject(io, readOnProc && writeOnProc);
236 
237  Info<< " Writing " << name << endl;
238  writeOk = meshObject.regIOobject::write(readOnProc && writeOnProc);
239  }
240 
241  return writeOk;
242 }
243 
244 
245 int main(int argc, char *argv[])
246 {
248  (
249  "Converts all IOobjects associated with a case into the format"
250  " specified in the controlDict"
251  );
254  (
255  "noConstant",
256  "Exclude the 'constant/' dir in the times list"
257  );
259  (
260  "enableFunctionEntries",
261  "Enable expansion of dictionary directives - #include, #codeStream etc"
262  );
263 
264  #include "addRegionOption.H"
265  #include "setRootCase.H"
266 
267  // enable noConstant by switching
268  if (!args.found("noConstant"))
269  {
270  args.setOption("constant");
271  }
272  else
273  {
274  args.unsetOption("constant");
275  Info<< "Excluding the constant directory." << nl << endl;
276  }
277 
278 
279  #include "createTime.H"
280  // Optional mesh (used to read Clouds)
282 
283  const bool enableEntries = args.found("enableFunctionEntries");
284  if (enableEntries)
285  {
286  Info<< "Allowing dictionary preprocessing ('#include', '#codeStream')."
287  << endl;
288  }
289 
290  const int oldFlag = entry::disableFunctionEntries;
291  if (!enableEntries)
292  {
293  // By default disable dictionary expansion for fields
295  }
296 
297  // Make sure we do not use the master-only reading since we read
298  // fields (different per processor) as dictionaries.
300 
301 
302  // Specified region or default region
303  #include "getRegionOption.H"
304  if (!polyMesh::regionName(regionName).empty())
305  {
306  Info<< "Using region " << regionName << nl << endl;
307  }
308 
309  const fileName meshDir
310  (
312  );
313 
314 
316 
317  forAll(timeDirs, timeI)
318  {
319  runTime.setTime(timeDirs[timeI], timeI);
320  Info<< "Time = " << runTime.timeName() << endl;
321 
322  // Convert all the standard mesh files
323  writeMeshObject<cellCompactIOList, cellIOList>
324  (
325  "cells",
326  meshDir,
327  runTime,
328  false // do not check typeName since varies between binary/ascii
329  );
330  writeMeshObject<labelIOList>("owner", meshDir, runTime);
331  writeMeshObject<labelIOList>("neighbour", meshDir, runTime);
332  writeMeshObject<faceCompactIOList, faceIOList>
333  (
334  "faces",
335  meshDir,
336  runTime,
337  false // do not check typeName since varies between binary/ascii
338  );
339  writeMeshObject<pointIOField>("points", meshDir, runTime);
340  // Write boundary in ascii. This is only needed for fileHandler to
341  // kick in. Should not give problems since always writing ascii.
342  writeZones("boundary", meshDir, runTime, IOstreamOption::UNCOMPRESSED);
343  writeMeshObject<labelIOList>("pointProcAddressing", meshDir, runTime);
344  writeMeshObject<labelIOList>("faceProcAddressing", meshDir, runTime);
345  writeMeshObject<labelIOList>("cellProcAddressing", meshDir, runTime);
346  writeMeshObject<labelIOList>
347  (
348  "boundaryProcAddressing",
349  meshDir,
350  runTime
351  );
352 
353  // foamyHexMesh vertices
354  writeMeshObject<pointIOField>
355  (
356  "internalDelaunayVertices",
358  runTime
359  );
360 
362  {
363  // Only do zones when converting from binary to ascii
364  // The other way gives problems since working on dictionary level.
365  const IOstreamOption::compressionType compress =
367  writeZones("cellZones", meshDir, runTime, compress);
368  writeZones("faceZones", meshDir, runTime, compress);
369  writeZones("pointZones", meshDir, runTime, compress);
370  }
371 
372  // Get list of objects from the database
373  IOobjectList objects
374  (
375  runTime,
376  runTime.timeName(),
378  );
379 
380  forAllConstIters(objects, iter)
381  {
382  const IOobject& io = *(iter.val());
383 
384  if
385  (
391 
397 
403 
409  )
410  {
411  Info<< " Reading " << io.headerClassName()
412  << " : " << io.name() << endl;
413 
415 
416  Info<< " Writing " << io.name() << endl;
417  fDict.regIOobject::write();
418  }
419  }
420 
421 
422  // Check for lagrangian
423  fileNameList lagrangianDirs
424  (
425  1,
426  fileHandler().filePath
427  (
428  runTime.timePath()
430  / cloud::prefix
431  )
432  );
433 
434  Pstream::combineReduce(lagrangianDirs, uniqueEqOp<fileName>());
435 
436  if (!lagrangianDirs.empty())
437  {
438  if (meshPtr)
439  {
440  meshPtr->readUpdate();
441  }
442  else
443  {
444  Info<< " Create polyMesh for time = "
445  << runTime.timeName() << endl;
446 
447  meshPtr.reset
448  (
449  new polyMesh
450  (
451  IOobject
452  (
454  runTime.timeName(),
455  runTime,
457  )
458  )
459  );
460  }
461 
462  const auto& mesh = meshPtr();
463 
464  fileNameList cloudDirs
465  (
467  (
468  lagrangianDirs[0],
470  )
471  );
472 
473  Pstream::combineReduce(cloudDirs, uniqueEqOp<fileName>());
474 
475  for (const auto& cloudDir : cloudDirs)
476  {
477  fileName dir(cloud::prefix/cloudDir);
478 
479  // Read with origProcId,origId fields
480  passiveParticleCloud parcels(mesh, cloudDir, true);
481 
482  const bool writeOnProc = parcels.size();
483 
484  parcels.writeObject
485  (
487  (
490  ),
491  writeOnProc
492  );
493 
494 
495  // Do local scan for valid cloud objects
496  wordList cloudFields
497  (
498  IOobjectList(runTime, runTime.timeName(), dir).sortedNames()
499  );
500 
501  // Combine with all other cloud objects
502  Pstream::combineReduce(cloudFields, uniqueEqOp<word>());
503 
504  for (const word& name : cloudFields)
505  {
506  // These ones already done by cloud itself
507  if
508  (
509  name == "positions"
510  || name == "coordinates"
511  || name == "origProcId"
512  || name == "origId"
513  )
514  {
515  continue;
516  }
517 
518  #undef doLocalCode
519  #define doLocalCode(Type) \
520  if \
521  ( \
522  writeOptionalMeshObject<Type> \
523  ( \
524  name, dir, runTime, writeOnProc \
525  ) \
526  ) \
527  { \
528  continue; \
529  }
530 
537 
540 
541  #undef doLocalCode
542 
543  Info<< " Failed converting " << name << endl;
544  }
545  }
546  }
547 
548  Info<< endl;
549  }
550 
552 
553  Info<< "End\n" << endl;
554 
555  return 0;
556 }
557 
558 
559 // ************************************************************************* //
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:367
static void addNote(const string &note)
Add extra notes for the usage information.
Definition: argList.C:462
A class for handling file names.
Definition: fileName.H:72
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:2241
A list of keyword definitions, which are a keyword followed by a number of values (eg...
Definition: dictionary.H:129
Write a mesh object in the format specified in controlDict.
fileName timePath() const
Return current time path = path/timeName.
Definition: Time.H:520
const word & regionName() const
The mesh region name or word::null if polyMesh::defaultRegion.
Definition: polyMesh.C:847
const word & name() const noexcept
Return the object name.
Definition: IOobjectI.H:195
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:50
"ascii" (normal default)
The meshObject is a concrete regIOobject to register MeshObject items.
Definition: MeshObject.H:90
engineTime & runTime
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
static void addBoolOption(const word &optName, const string &usage="", bool advanced=false)
Add a bool option to validOptions with usage information.
Definition: argList.C:374
refPtr< fileOperation > fileHandler(std::nullptr_t)
Delete current file handler - forwards to fileOperation::handler()
Required Classes.
A simple container for options an IOstream can normally have.
Ignore writing from objectRegistry::writeObject()
#define doLocalCode(FieldType, Variable)
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:69
A Cloud of passive particles.
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry (const access) with the given keyword.
Definition: dictionaryI.H:104
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
scalar y
Required Classes.
static int disableFunctionEntries
Enable or disable use of function entries and variable expansions.
Definition: entry.H:139
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
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:406
static void combineReduce(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...
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:902
bool typeHeaderOk(const bool checkType=true, const bool search=true, const bool verbose=true)
Read header (respects is_globalIOobject trait) and check its info.
Read field as dictionary (without mesh).
IOstreamOption::compressionType writeCompression() const noexcept
Get the write stream compression.
Definition: TimeI.H:143
static word timeName(const scalar t, const int precision=precision_)
Return a time name for the given scalar time value formatted with the given precision.
Definition: Time.C:714
defineTemplateTypeNameAndDebug(faScalarMatrix, 0)
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:351
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:260
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:213
Foam::word regionName(args.getOrDefault< word >("region", Foam::polyMesh::defaultRegion))
bool setOption(const word &optName, const string &param="")
Set option directly (use with caution)
Definition: argList.C:2213
fileName meshDir() const
Return the local mesh directory (dbDir()/meshSubDir)
Definition: polyMesh.C:841
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
regIOobject is an abstract class derived from IOobject to handle automatic object registration with t...
Definition: regIOobject.H:68
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:75
entry * set(entry *entryPtr)
Assign a new entry, overwriting any existing entry.
Definition: dictionary.C:765
IOstreamOption::streamFormat writeFormat() const noexcept
Get write stream format.
Definition: TimeI.H:123
List< label > labelList
A List of labels.
Definition: List.H:62
IOobject io("surfaceFilmProperties", mesh.time().constant(), mesh, IOobject::READ_IF_PRESENT, IOobject::NO_WRITE, IOobject::NO_REGISTER)
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:963
Foam::argList args(argc, argv)
bool returnReduceOr(const bool value, const label comm=UPstream::worldComm)
Perform logical (or) MPI Allreduce on a copy. Uses UPstream::reduceOr.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
List< bool > boolList
A List of bools.
Definition: List.H:60
A primitive field of type <T> with automated input and output.
Do not request registration (bool: false)
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
A keyword and a list of tokens is an &#39;entry&#39;.
Definition: entry.H:63
bool isHeaderClass() const
Check if headerClassName() equals Type::typeName.
Definition: IOobjectI.H:258
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:79