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 | Copyright (C) 2011-2015 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 
53 #include "writeMeshObject.H"
54 #include "fieldDictionary.H"
55 
56 using namespace Foam;
57 
58 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
59 
60 namespace Foam
61 {
63 }
64 
65 
66 // Hack to do zones which have Lists in them. See above.
67 bool writeZones(const word& name, const fileName& meshDir, Time& runTime)
68 {
69  IOobject io
70  (
71  name,
72  runTime.timeName(),
73  meshDir,
74  runTime,
77  false
78  );
79 
80  bool writeOk = false;
81 
82  if (io.headerOk())
83  {
84  Info<< " Reading " << io.headerClassName()
85  << " : " << name << endl;
86 
87  // Switch off type checking (for reading e.g. faceZones as
88  // generic list of dictionaries).
89  const word oldTypeName = IOPtrList<entry>::typeName;
91 
93 
94  forAll(meshObject, i)
95  {
96  if (meshObject[i].isDict())
97  {
98  dictionary& d = meshObject[i].dict();
99 
100  if (d.found("faceLabels"))
101  {
102  d.set("faceLabels", labelList(d.lookup("faceLabels")));
103  }
104 
105  if (d.found("flipMap"))
106  {
107  d.set("flipMap", boolList(d.lookup("flipMap")));
108  }
109 
110  if (d.found("cellLabels"))
111  {
112  d.set("cellLabels", labelList(d.lookup("cellLabels")));
113  }
114 
115  if (d.found("pointLabels"))
116  {
117  d.set("pointLabels", labelList(d.lookup("pointLabels")));
118  }
119  }
120  }
121 
122  const_cast<word&>(IOPtrList<entry>::typeName) = oldTypeName;
123  // Fake type back to what was in field
124  const_cast<word&>(meshObject.type()) = io.headerClassName();
125 
126  Info<< " Writing " << name << endl;
127 
128  // Force writing as ascii
129  writeOk = meshObject.regIOobject::writeObject
130  (
133  runTime.writeCompression()
134  );
135  }
136 
137  return writeOk;
138 }
139 
140 
141 
142 
143 int main(int argc, char *argv[])
144 {
147  (
148  "noConstant",
149  "exclude the 'constant/' dir in the times list"
150  );
151 
152  #include "addRegionOption.H"
153  #include "setRootCase.H"
154 
155  // enable noConstant by switching
156  if (!args.optionFound("noConstant"))
157  {
158  args.setOption("constant", "");
159  }
160  else
161  {
162  args.unsetOption("constant");
163  Info<< "Excluding the constant directory." << nl << endl;
164  }
165 
166 
167  #include "createTime.H"
168 
169 
170  // Make sure we do not use the master-only reading since we read
171  // fields (different per processor) as dictionaries.
173 
174 
175  fileName meshDir = polyMesh::meshSubDir;
176  fileName regionPrefix = "";
177  word regionName = polyMesh::defaultRegion;
178  if (args.optionReadIfPresent("region", regionName))
179  {
180  Info<< "Using region " << regionName << nl << endl;
181  regionPrefix = regionName;
182  meshDir = regionName/polyMesh::meshSubDir;
183  }
184 
186 
187  forAll(timeDirs, timeI)
188  {
189  runTime.setTime(timeDirs[timeI], timeI);
190  Info<< "Time = " << runTime.timeName() << endl;
191 
192  // Convert all the standard mesh files
193  writeMeshObject<cellCompactIOList>("cells", meshDir, runTime);
194  writeMeshObject<labelIOList>("owner", meshDir, runTime);
195  writeMeshObject<labelIOList>("neighbour", meshDir, runTime);
196  writeMeshObject<faceCompactIOList>("faces", meshDir, runTime);
197  writeMeshObject<pointIOField>("points", meshDir, runTime);
198  writeMeshObject<labelIOList>("pointProcAddressing", meshDir, runTime);
199  writeMeshObject<labelIOList>("faceProcAddressing", meshDir, runTime);
200  writeMeshObject<labelIOList>("cellProcAddressing", meshDir, runTime);
201  writeMeshObject<labelIOList>
202  (
203  "boundaryProcAddressing",
204  meshDir,
205  runTime
206  );
207 
208  // foamyHexMesh vertices
209  writeMeshObject<pointIOField>
210  (
211  "internalDelaunayVertices",
212  regionPrefix,
213  runTime
214  );
215 
216  if (runTime.writeFormat() == IOstream::ASCII)
217  {
218  // Only do zones when converting from binary to ascii
219  // The other way gives problems since working on dictionary level.
220  writeZones("cellZones", meshDir, runTime);
221  writeZones("faceZones", meshDir, runTime);
222  writeZones("pointZones", meshDir, runTime);
223  }
224 
225  // Get list of objects from the database
226  IOobjectList objects(runTime, runTime.timeName(), regionPrefix);
227 
228  forAllConstIter(IOobjectList, objects, iter)
229  {
230  const word& headerClassName = iter()->headerClassName();
231 
232  if
233  (
234  headerClassName == volScalarField::typeName
235  || headerClassName == volVectorField::typeName
236  || headerClassName == volSphericalTensorField::typeName
237  || headerClassName == volSymmTensorField::typeName
238  || headerClassName == volTensorField::typeName
239 
240  || headerClassName == surfaceScalarField::typeName
241  || headerClassName == surfaceVectorField::typeName
242  || headerClassName == surfaceSphericalTensorField::typeName
243  || headerClassName == surfaceSymmTensorField::typeName
244  || headerClassName == surfaceTensorField::typeName
245 
246  || headerClassName == pointScalarField::typeName
247  || headerClassName == pointVectorField::typeName
248  || headerClassName == pointSphericalTensorField::typeName
249  || headerClassName == pointSymmTensorField::typeName
250  || headerClassName == pointTensorField::typeName
251  )
252  {
253  Info<< " Reading " << headerClassName
254  << " : " << iter()->name() << endl;
255 
256  fieldDictionary fDict
257  (
258  *iter(),
259  headerClassName
260  );
261 
262  Info<< " Writing " << iter()->name() << endl;
263  fDict.regIOobject::write();
264  }
265  }
266 
267  Info<< endl;
268  }
269 
270  Info<< "End\n" << endl;
271 
272  return 0;
273 }
274 
275 
276 // ************************************************************************* //
Foam::surfaceFields.
bool setOption(const word &opt, const string &param="")
Set option directly (use with caution)
Definition: argList.C:909
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:428
A class for handling file names.
Definition: fileName.H:69
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:50
bool unsetOption(const word &opt)
Unset option directly (use with caution)
Definition: argList.C:984
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 word meshSubDir
Return the mesh sub-directory name (usually "polyMesh")
Definition: polyMesh.H:309
bool optionReadIfPresent(const word &opt, T &) const
Read a value from the named option if present.
Definition: argListI.H:198
static word defaultRegion
Return the default region name.
Definition: polyMesh.H:306
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:253
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:715
Class to control time during OpenFOAM simulations that is also the top-level objectRegistry.
Definition: Time.H:68
List< bool > boolList
Bool container classes.
Definition: boolList.H:50
static fileCheckTypes fileModificationChecking
Definition: regIOobject.H:128
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
A class for handling words, derived from string.
Definition: word.H:59
static const word null
An empty word.
Definition: word.H:77
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:306
List< label > labelList
A List of labels.
Definition: labelList.H:56
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:262
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:253
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:864
messageStream Info
static void addBoolOption(const word &opt, const string &usage="")
Add to a bool option to validOptions with usage information.
Definition: argList.C:83
Foam::argList args(argc, argv)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
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:451