initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / miscellaneous / foamFormatConvert / foamFormatConvert.C
blob0140ba4e232affb451fc499b1ba2048fa7361ab0
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     foamFormatConvert
28 Description
29     Converts all IOobjects associated with a case into the format specified
30     in the controlDict.
32     Mainly used to convert binary mesh/field files to ASCII.
34 \*---------------------------------------------------------------------------*/
36 #include "argList.H"
37 #include "timeSelector.H"
38 #include "Time.H"
39 #include "volFields.H"
40 #include "surfaceFields.H"
41 #include "pointFields.H"
42 #include "cellIOList.H"
43 #include "IOobjectList.H"
44 #include "IOPtrList.H"
46 #include "writeMeshObject.H"
47 #include "fieldDictionary.H"
49 using namespace Foam;
51 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
53 namespace Foam
55     defineTemplateTypeNameAndDebug(IOPtrList<entry>, 0);
59 // Main program:
61 int main(int argc, char *argv[])
63     timeSelector::addOptions();
64 #   include "setRootCase.H"
65 #   include "createTime.H"
66     Foam::instantList timeDirs = Foam::timeSelector::select0(runTime, args);
68     forAll (timeDirs, timeI)
69     {
70         runTime.setTime(timeDirs[timeI], timeI);
71         Info<< "Time = " << runTime.timeName() << endl;
73         // Convert all the standard mesh files
74         writeMeshObject<cellIOList>("cells", runTime);
75         writeMeshObject<labelIOList>("owner", runTime);
76         writeMeshObject<labelIOList>("neighbour", runTime);
77         writeMeshObject<faceIOList>("faces", runTime);
78         writeMeshObject<pointIOField>("points", runTime);
79         writeMeshObject<IOPtrList<entry> >("cellZones", runTime);
80         writeMeshObject<IOPtrList<entry> >("faceZones", runTime);
81         writeMeshObject<IOPtrList<entry> >("pointZones", runTime);
83         // Get list of objects from the database
84         IOobjectList objects(runTime, runTime.timeName());
86         forAllConstIter(IOobjectList, objects, iter)
87         {
88             const word& headerClassName = iter()->headerClassName();
90             if
91             (
92                 headerClassName == volScalarField::typeName
93              || headerClassName == volVectorField::typeName
94              || headerClassName == volSphericalTensorField::typeName
95              || headerClassName == volSymmTensorField::typeName
96              || headerClassName == volTensorField::typeName
98              || headerClassName == surfaceScalarField::typeName
99              || headerClassName == surfaceVectorField::typeName
100              || headerClassName == surfaceSphericalTensorField::typeName
101              || headerClassName == surfaceSymmTensorField::typeName
102              || headerClassName == surfaceTensorField::typeName
104              || headerClassName == pointScalarField::typeName
105              || headerClassName == pointVectorField::typeName
106              || headerClassName == pointSphericalTensorField::typeName
107              || headerClassName == pointSymmTensorField::typeName
108              || headerClassName == pointTensorField::typeName
109             )
110             {
111                 Info<< "        Reading " << headerClassName
112                     << " : " << iter()->name() << endl;
114                 fieldDictionary fDict
115                 (
116                     *iter(),
117                     headerClassName
118                 );
120                 Info<< "        Writing " << iter()->name() << endl;
121                 fDict.regIOobject::write();
122             }
123         }
125         Info<< endl;
126     }
128     Info<< "End\n" << endl;
130     return 0;
134 // ************************************************************************* //