BUG: UListIO: byteSize overflowing on really big faceLists
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / dataConversion / foamToEnsight / foamToEnsight.C
blobe9a942f38a24b47796a74d9e9db4d94a304cf314
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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.
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.
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/>.
24 Description
25     Translates OpenFOAM data to EnSight format.
27     An Ensight part is created for the internalMesh and for each patch.
29 Usage
30     - foamToEnsight [OPTION] \n
31     Translates OpenFOAM data to Ensight format
33     \param -ascii \n
34     Write Ensight data in ASCII format instead of "C Binary"
36     \param -patches patchList \n
37     Specify particular patches to write.
38     Specifying an empty list suppresses writing the internalMesh.
40     \param -noPatches \n
41     Suppress writing any patches.
43     \param -faceZones zoneList \n
44     Specify faceZones to write, with wildcards
46 Note
47     Parallel support for cloud data is not supported
48     - writes to \a EnSight directory to avoid collisions with foamToEnsightParts
50 \*---------------------------------------------------------------------------*/
52 #include "argList.H"
53 #include "timeSelector.H"
54 #include "IOobjectList.H"
55 #include "IOmanip.H"
56 #include "OFstream.H"
58 #include "volFields.H"
60 #include "labelIOField.H"
61 #include "scalarIOField.H"
62 #include "tensorIOField.H"
64 #include "ensightMesh.H"
65 #include "ensightField.H"
67 #include "ensightParticlePositions.H"
68 #include "ensightCloudField.H"
70 #include "fvc.H"
72 using namespace Foam;
74 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
76 bool inFileNameList
78     const fileNameList& nameList,
79     const word& name
82     forAll(nameList, i)
83     {
84         if (nameList[i] == name)
85         {
86             return true;
87         }
88     }
90     return false;
94 // Main program:
96 int main(int argc, char *argv[])
98     timeSelector::addOptions();
99 #   include "addRegionOption.H"
101     argList::addBoolOption
102     (
103         "ascii",
104         "write in ASCII format instead of 'C Binary'"
105     );
106     argList::addBoolOption
107     (
108         "nodeValues",
109         "write values in nodes"
110     );
111     argList::addBoolOption
112     (
113         "noPatches",
114         "suppress writing any patches"
115     );
116     argList::addOption
117     (
118         "patches",
119         "wordReList",
120         "specify particular patches to write - eg '(outlet \"inlet.*\")'. "
121         "An empty list suppresses writing the internalMesh."
122     );
123     argList::addOption
124     (
125         "faceZones",
126         "wordReList",
127         "specify faceZones to write - eg '( slice \"mfp-.*\" )'."
128     );
130 #   include "setRootCase.H"
132     // Check options
133     const bool binary = !args.optionFound("ascii");
134     const bool nodeValues = args.optionFound("nodeValues");
136 #   include "createTime.H"
138     instantList Times = timeSelector::select0(runTime, args);
140 #   include "createNamedMesh.H"
142     // Mesh instance (region0 gets filtered out)
143     fileName regionPrefix = "";
145     if (regionName != polyMesh::defaultRegion)
146     {
147         regionPrefix = regionName;
148     }
150     const label nVolFieldTypes = 5;
151     const word volFieldTypes[] =
152     {
153         volScalarField::typeName,
154         volVectorField::typeName,
155         volSphericalTensorField::typeName,
156         volSymmTensorField::typeName,
157         volTensorField::typeName
158     };
160     // Path to EnSight folder at case level only
161     // - For parallel cases, data only written from master
162     fileName ensightDir = args.rootPath()/args.globalCaseName()/"EnSight";
164     if (Pstream::master())
165     {
166         if (isDir(ensightDir))
167         {
168             rmDir(ensightDir);
169         }
171         mkDir(ensightDir);
172     }
174     // Start of case file header output
175     // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~~
177     const word prepend = args.globalCaseName() + '.';
179     OFstream *ensightCaseFilePtr = NULL;
180     if (Pstream::master())
181     {
182         fileName caseFileName = prepend + "case";
183         Info<< nl << "write case: " << caseFileName.c_str() << endl;
185         // the case file is always ASCII
186         ensightCaseFilePtr = new OFstream
187         (
188             ensightDir/caseFileName,
189             IOstream::ASCII
190         );
192         *ensightCaseFilePtr
193             << "FORMAT" << nl
194             << "type: ensight gold" << nl << nl;
195     }
197     OFstream& ensightCaseFile = *ensightCaseFilePtr;
199     // Construct the EnSight mesh
200     const bool selectedPatches = args.optionFound("patches");
201     wordReList patchPatterns;
202     if (selectedPatches)
203     {
204         patchPatterns = wordReList(args.optionLookup("patches")());
205     }
206     const bool selectedZones = args.optionFound("faceZones");
207     wordReList zonePatterns;
208     if (selectedZones)
209     {
210         zonePatterns = wordReList(args.optionLookup("faceZones")());
211     }
213     ensightMesh eMesh
214     (
215         mesh,
216         args.optionFound("noPatches"),
217         selectedPatches,
218         patchPatterns,
219         selectedZones,
220         zonePatterns,
221         binary
222     );
224     // Set Time to the last time before looking for the lagrangian objects
225     runTime.setTime(Times.last(), Times.size()-1);
227     IOobjectList objects(mesh, runTime.timeName());
229 #   include "checkMeshMoving.H"
231     wordHashSet allCloudNames;
232     if (Pstream::master())
233     {
234         word geomFileName = prepend + "000";
236         // test pre check variable if there is a moving mesh
237         if (meshMoving)
238         {
239             geomFileName = prepend + "***";
240         }
242         ensightCaseFile
243             << "GEOMETRY" << nl
244             << "model:        1     "
245             << (geomFileName + ".mesh").c_str() << nl;
246     }
248     // Identify if lagrangian data exists at each time, and add clouds
249     // to the 'allCloudNames' hash set
250     forAll(Times, timeI)
251     {
252         runTime.setTime(Times[timeI], timeI);
254         fileNameList cloudDirs = readDir
255         (
256             runTime.timePath()/regionPrefix/cloud::prefix,
257             fileName::DIRECTORY
258         );
260         forAll(cloudDirs, cloudI)
261         {
262             IOobjectList cloudObjs
263             (
264                 mesh,
265                 runTime.timeName(),
266                 cloud::prefix/cloudDirs[cloudI]
267             );
269             IOobject* positionsPtr = cloudObjs.lookup("positions");
271             if (positionsPtr)
272             {
273                 allCloudNames.insert(cloudDirs[cloudI]);
274             }
275         }
276     }
278     HashTable<HashTable<word> > allCloudFields;
279     forAllConstIter(wordHashSet, allCloudNames, cloudIter)
280     {
281         // Add the name of the cloud(s) to the case file header
282         if (Pstream::master())
283         {
284             ensightCaseFile
285             <<  (
286                     "measured:     1     "
287                   + prepend
288                   + "***."
289                   + cloudIter.key()
290                 ).c_str()
291             << nl;
292         }
294         // Create a new hash table for each cloud
295         allCloudFields.insert(cloudIter.key(), HashTable<word>());
297         // Identify the new cloud in the hash table
298         HashTable<HashTable<word> >::iterator newCloudIter =
299             allCloudFields.find(cloudIter.key());
301         // Loop over all times to build list of fields and field types
302         // for each cloud
303         forAll(Times, timeI)
304         {
305             runTime.setTime(Times[timeI], timeI);
307             IOobjectList cloudObjs
308             (
309                 mesh,
310                 runTime.timeName(),
311                 cloud::prefix/cloudIter.key()
312             );
314             forAllConstIter(IOobjectList, cloudObjs, fieldIter)
315             {
316                 const IOobject obj = *fieldIter();
318                 if (obj.name() != "positions")
319                 {
320                     // Add field and field type
321                     newCloudIter().insert
322                     (
323                         obj.name(),
324                         obj.headerClassName()
325                     );
326                 }
327             }
328         }
329     }
331     label nTimeSteps = 0;
332     forAll(Times, timeIndex)
333     {
334         nTimeSteps++;
335         runTime.setTime(Times[timeIndex], timeIndex);
337         word timeName = itoa(timeIndex);
338         word timeFile = prepend + timeName;
340         Info<< "Translating time = " << runTime.timeName() << nl;
342         polyMesh::readUpdateState meshState = mesh.readUpdate();
344         if (meshState != polyMesh::UNCHANGED)
345         {
346             eMesh.correct();
347         }
349         if (timeIndex == 0 || (meshState != polyMesh::UNCHANGED))
350         {
351             eMesh.write
352             (
353                 ensightDir,
354                 prepend,
355                 timeIndex,
356                 ensightCaseFile
357             );
358         }
361         // Start of field data output
362         // ~~~~~~~~~~~~~~~~~~~~~~~~~~
364         if (timeIndex == 0 && Pstream::master())
365         {
366             ensightCaseFile<< nl << "VARIABLE" << nl;
367         }
370         // Cell field data output
371         // ~~~~~~~~~~~~~~~~~~~~~~
373         for (label i=0; i<nVolFieldTypes; i++)
374         {
375             wordList fieldNames = objects.names(volFieldTypes[i]);
377             forAll(fieldNames, j)
378             {
379                 const word& fieldName = fieldNames[j];
381 #               include "checkData.H"
383                 if (!variableGood)
384                 {
385                     continue;
386                 }
388                 IOobject fieldObject
389                 (
390                     fieldName,
391                     mesh.time().timeName(),
392                     mesh,
393                     IOobject::MUST_READ,
394                     IOobject::NO_WRITE
395                 );
397                 if (volFieldTypes[i] == volScalarField::typeName)
398                 {
399                     ensightField<scalar>
400                     (
401                         fieldObject,
402                         eMesh,
403                         ensightDir,
404                         prepend,
405                         timeIndex,
406                         binary,
407                         nodeValues,
408                         ensightCaseFile
409                     );
410                 }
411                 else if (volFieldTypes[i] == volVectorField::typeName)
412                 {
413                     ensightField<vector>
414                     (
415                         fieldObject,
416                         eMesh,
417                         ensightDir,
418                         prepend,
419                         timeIndex,
420                         binary,
421                         nodeValues,
422                         ensightCaseFile
423                     );
424                 }
425                 else if (volFieldTypes[i] == volSphericalTensorField::typeName)
426                 {
427                     ensightField<sphericalTensor>
428                     (
429                         fieldObject,
430                         eMesh,
431                         ensightDir,
432                         prepend,
433                         timeIndex,
434                         binary,
435                         nodeValues,
436                         ensightCaseFile
437                     );
438                 }
439                 else if (volFieldTypes[i] == volSymmTensorField::typeName)
440                 {
441                     ensightField<symmTensor>
442                     (
443                         fieldObject,
444                         eMesh,
445                         ensightDir,
446                         prepend,
447                         timeIndex,
448                         binary,
449                         nodeValues,
450                         ensightCaseFile
451                     );
452                 }
453                 else if (volFieldTypes[i] == volTensorField::typeName)
454                 {
455                     ensightField<tensor>
456                     (
457                         fieldObject,
458                         eMesh,
459                         ensightDir,
460                         prepend,
461                         timeIndex,
462                         binary,
463                         nodeValues,
464                         ensightCaseFile
465                     );
466                 }
467             }
468         }
471         // Cloud field data output
472         // ~~~~~~~~~~~~~~~~~~~~~~~
474         forAllConstIter(HashTable<HashTable<word> >, allCloudFields, cloudIter)
475         {
476             const word& cloudName = cloudIter.key();
478             fileNameList currentCloudDirs = readDir
479             (
480                 runTime.timePath()/regionPrefix/cloud::prefix,
481                 fileName::DIRECTORY
482             );
484             bool cloudExists = inFileNameList(currentCloudDirs, cloudName);
485             ensightParticlePositions
486             (
487                 mesh,
488                 ensightDir,
489                 timeFile,
490                 cloudName,
491                 cloudExists
492             );
494             forAllConstIter(HashTable<word>, cloudIter(), fieldIter)
495             {
496                 const word& fieldName = fieldIter.key();
497                 const word& fieldType = fieldIter();
499                 IOobject fieldObject
500                 (
501                     fieldName,
502                     mesh.time().timeName(),
503                     cloud::prefix/cloudName,
504                     mesh,
505                     IOobject::MUST_READ
506                 );
508                 bool fieldExists = fieldObject.headerOk();
509                 if (fieldType == scalarIOField::typeName)
510                 {
511                     ensightCloudField<scalar>
512                     (
513                         fieldObject,
514                         ensightDir,
515                         prepend,
516                         timeIndex,
517                         cloudName,
518                         ensightCaseFile,
519                         fieldExists
520                     );
521                 }
522                 else if (fieldType == vectorIOField::typeName)
523                 {
524                     ensightCloudField<vector>
525                     (
526                         fieldObject,
527                         ensightDir,
528                         prepend,
529                         timeIndex,
530                         cloudName,
531                         ensightCaseFile,
532                         fieldExists
533                     );
534                 }
535                 else
536                 {
537                     Info<< "Unable to convert field type " << fieldType
538                         << " for field " << fieldName << endl;
539                 }
540             }
541         }
542     }
544 #   include "ensightCaseTail.H"
546     if (Pstream::master())
547     {
548         delete ensightCaseFilePtr;
549     }
551     Info<< "End\n" << endl;
553     return 0;
557 // ************************************************************************* //