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-2018 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.typeHeaderOk<cellZoneMesh>(false))
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 valid
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.typeHeaderOk<IOField<label>>(false);
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, valid && haveFile);
214 
215  Info<< " Writing " << name << endl;
216  writeOk = meshObject.regIOobject::write(valid && 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.timeName() << 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  writeMeshObject<labelIOList>
297  (
298  "boundaryProcAddressing",
299  meshDir,
300  runTime
301  );
302 
303  // foamyHexMesh vertices
304  writeMeshObject<pointIOField>
305  (
306  "internalDelaunayVertices",
307  regionPrefix,
308  runTime
309  );
310 
311  if (runTime.writeFormat() == IOstream::ASCII)
312  {
313  // Only do zones when converting from binary to ascii
314  // The other way gives problems since working on dictionary level.
315  writeZones("cellZones", meshDir, runTime);
316  writeZones("faceZones", meshDir, runTime);
317  writeZones("pointZones", meshDir, runTime);
318  }
319 
320  // Get list of objects from the database
321  IOobjectList objects(runTime, runTime.timeName(), regionPrefix);
322 
323  forAllConstIter(IOobjectList, objects, iter)
324  {
325  const word& headerClassName = iter()->headerClassName();
326 
327  if
328  (
329  headerClassName == volScalarField::typeName
330  || headerClassName == volVectorField::typeName
331  || headerClassName == volSphericalTensorField::typeName
332  || headerClassName == volSymmTensorField::typeName
333  || headerClassName == volTensorField::typeName
334 
335  || headerClassName == surfaceScalarField::typeName
336  || headerClassName == surfaceVectorField::typeName
337  || headerClassName == surfaceSphericalTensorField::typeName
338  || headerClassName == surfaceSymmTensorField::typeName
339  || headerClassName == surfaceTensorField::typeName
340 
341  || headerClassName == pointScalarField::typeName
342  || headerClassName == pointVectorField::typeName
343  || headerClassName == pointSphericalTensorField::typeName
344  || headerClassName == pointSymmTensorField::typeName
345  || headerClassName == pointTensorField::typeName
346 
347  || headerClassName == volScalarField::Internal::typeName
348  || headerClassName == volVectorField::Internal::typeName
349  || headerClassName == volSphericalTensorField::Internal::typeName
350  || headerClassName == volSymmTensorField::Internal::typeName
351  || headerClassName == volTensorField::Internal::typeName
352  )
353  {
354  Info<< " Reading " << headerClassName
355  << " : " << iter()->name() << endl;
356 
357  fieldDictionary fDict(*iter(), headerClassName);
358 
359  Info<< " Writing " << iter()->name() << endl;
360  fDict.regIOobject::write();
361  }
362  }
363 
364 
365 
366  // Check for lagrangian
367  stringList lagrangianDirs
368  (
369  1,
370  fileHandler().filePath
371  (
372  runTime.timePath()
373  / regionPrefix
374  / cloud::prefix
375  )
376  );
377 
378  combineReduce(lagrangianDirs, uniqueEqOp());
379 
380  if (!lagrangianDirs.empty())
381  {
382  if (meshPtr.valid())
383  {
384  meshPtr().readUpdate();
385  }
386  else
387  {
388  Info<< " Create polyMesh for time = "
389  << runTime.timeName() << endl;
390 
391  meshPtr.reset
392  (
393  new polyMesh
394  (
395  IOobject
396  (
398  runTime.timeName(),
399  runTime,
401  )
402  )
403  );
404  }
405 
406  stringList cloudDirs
407  (
409  (
410  lagrangianDirs[0],
412  )
413  );
414 
415  combineReduce(cloudDirs, uniqueEqOp());
416 
417  forAll(cloudDirs, i)
418  {
419  fileName dir(cloud::prefix/cloudDirs[i]);
420 
421  Cloud<passiveParticle> parcels(meshPtr(), cloudDirs[i], false);
422 
423  parcels.writeObject
424  (
425  runTime.writeFormat(),
427  runTime.writeCompression(),
428  parcels.size()
429  );
430 
431 
432  // Do local scan for valid cloud objects
433  IOobjectList sprayObjs(runTime, runTime.timeName(), dir);
434 
435  // Combine with all other cloud objects
436  stringList sprayFields(sprayObjs.sortedToc());
437  combineReduce(sprayFields, uniqueEqOp());
438 
439  forAll(sprayFields, fieldi)
440  {
441  const word& name = sprayFields[fieldi];
442 
443  // Note: try the various field types. Make sure to
444  // exit once successful conversion to avoid re-read
445  // converted file.
446 
447  if
448  (
449  name == "positions"
450  || name == "origProcId"
451  || name == "origId"
452  )
453  {
454  continue;
455  }
456 
457  bool writeOk = writeOptionalMeshObject<labelIOField>
458  (
459  name,
460  dir,
461  runTime,
462  parcels.size() > 0
463  );
464  if (writeOk) continue;
465 
466  writeOk = writeOptionalMeshObject<scalarIOField>
467  (
468  name,
469  dir,
470  runTime,
471  parcels.size() > 0
472  );
473  if (writeOk) continue;
474 
475  writeOk = writeOptionalMeshObject<vectorIOField>
476  (
477  name,
478  dir,
479  runTime,
480  parcels.size() > 0
481  );
482  if (writeOk) continue;
483 
484  writeOk = writeOptionalMeshObject<sphericalTensorIOField>
485  (
486  name,
487  dir,
488  runTime,
489  parcels.size() > 0
490  );
491  if (writeOk) continue;
492 
493  writeOk = writeOptionalMeshObject<symmTensorIOField>
494  (
495  name,
496  dir,
497  runTime,
498  parcels.size() > 0
499  );
500  if (writeOk) continue;
501 
502  writeOk = writeOptionalMeshObject<tensorIOField>
503  (
504  name,
505  dir,
506  runTime,
507  parcels.size() > 0
508  );
509  if (writeOk) continue;
510 
511  writeOk = writeOptionalMeshObject<labelFieldIOField>
512  (
513  name,
514  dir,
515  runTime,
516  parcels.size() > 0
517  );
518  if (writeOk) continue;
519 
520  writeOk = writeOptionalMeshObject<vectorFieldIOField>
521  (
522  name,
523  dir,
524  runTime,
525  parcels.size() > 0
526  );
527 
528  if (!writeOk)
529  {
530  Info<< " Failed converting " << name << endl;
531  }
532  }
533  }
534  }
535 
536  Info<< endl;
537  }
538 
539  Info<< "End\n" << endl;
540 
541  return 0;
542 }
543 
544 
545 // ************************************************************************* //
Foam::surfaceFields.
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:431
bool setOption(const word &opt, const string &param="")
Set option directly (use with caution)
Definition: argList.C:1064
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
void combineReduce(const List< UPstream::commsStruct > &comms, T &Value, const CombineOp &cop, const int tag, const label comm)
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
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
A class for handling file names.
Definition: fileName.H:69
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:1140
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:137
Write a mesh object in the format specified in controlDict.
static const char *const typeName
Definition: Field.H:94
Foam::word regionName
static fvMesh * meshPtr
Definition: globalFoam.H:52
static word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:312
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
engineTime & runTime
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:309
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
label size() const
Return number of elements in table.
Definition: HashTableI.H:65
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
virtual readUpdateState readUpdate()
Update the mesh based on the mesh files saved in time.
Definition: fvMesh.C:494
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.
bool valid() const
Return true if the autoPtr valid (ie, the pointer is set)
Definition: autoPtrI.H:83
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
defineTemplateTypeNameAndDebug(IOPtrList< ensightPart >, 0)
Read field as dictionary (without mesh).
static const char nl
Definition: Ostream.H:265
static fileCheckTypes fileModificationChecking
Type of file modification checking.
Definition: IOobject.H:214
static const word prefix
The prefix to local: lagrangian.
Definition: cloud.H:71
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:257
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
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:206
void set(entry *)
Assign a new entry, overwrite any existing entry.
Definition: dictionary.C:941
fileNameList readDir(const fileName &, const fileName::Type=fileName::FILE, const bool filtergz=true, const bool followLink=true)
Read a directory and return the entries as a string list.
Definition: POSIX.C:635
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:74
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:110
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
A primitive field of type <T> with automated input and output.
Definition: IOField.H:50
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:576