Initial commit for version 2.0.x patch release
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / dataConversion / foamDataToFluent / foamDataToFluent.C
blobb1ff214a772f84229d14e0ece7501be2c88f838c
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2004-2010 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
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 Fluent format.
27 \*---------------------------------------------------------------------------*/
29 #include "fvCFD.H"
30 #include "writeFluentFields.H"
31 #include "OFstream.H"
32 #include "IOobjectList.H"
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35 // Main program:
37 int main(int argc, char *argv[])
39     argList::noParallel();
40     timeSelector::addOptions(false);   // no constant
42 #   include "setRootCase.H"
43 #   include "createTime.H"
45     instantList timeDirs = timeSelector::select0(runTime, args);
47 #   include "createMesh.H"
49     // make a directory called proInterface in the case
50     mkDir(runTime.rootPath()/runTime.caseName()/"fluentInterface");
52     forAll(timeDirs, timeI)
53     {
54         runTime.setTime(timeDirs[timeI], timeI);
56         Info<< "Time = " << runTime.timeName() << endl;
58         if (mesh.readUpdate())
59         {
60             Info<< "    Read new mesh" << endl;
61         }
63         // make a directory called proInterface in the case
64         mkDir(runTime.rootPath()/runTime.caseName()/"fluentInterface");
66         // open a file for the mesh
67         OFstream fluentDataFile
68         (
69             runTime.rootPath()/
70             runTime.caseName()/
71             "fluentInterface"/
72             runTime.caseName() + runTime.timeName() + ".dat"
73         );
75         fluentDataFile
76             << "(0 \"FOAM to Fluent data File\")" << endl << endl;
78         // Writing number of faces
79         label nFaces = mesh.nFaces();
81         forAll(mesh.boundary(), patchI)
82         {
83             nFaces += mesh.boundary()[patchI].size();
84         }
86         fluentDataFile
87             << "(33 (" << mesh.nCells() << " " << nFaces << " "
88             << mesh.nPoints() << "))" << endl;
90         IOdictionary foamDataToFluentDict
91         (
92             IOobject
93             (
94                 "foamDataToFluentDict",
95                 runTime.system(),
96                 mesh,
97                 IOobject::MUST_READ_IF_MODIFIED,
98                 IOobject::NO_WRITE
99             )
100         );
103         // Search for list of objects for this time
104         IOobjectList objects(mesh, runTime.timeName());
107         // Converting volScalarField
108         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
110         // Search list of objects for volScalarFields
111         IOobjectList scalarFields(objects.lookupClass("volScalarField"));
113         forAllIter(IOobjectList, scalarFields, iter)
114         {
115             // Read field
116             volScalarField field(*iter(), mesh);
118             // lookup field from dictionary and convert field
119             label unitNumber;
120             if
121             (
122                 foamDataToFluentDict.readIfPresent(field.name(), unitNumber)
123              && unitNumber > 0
124             )
125             {
126                 Info<< "    Converting field " << field.name() << endl;
127                 writeFluentField(field, unitNumber, fluentDataFile);
128             }
129         }
132         // Converting volVectorField
133         // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
135         // Search list of objects for volVectorFields
136         IOobjectList vectorFields(objects.lookupClass("volVectorField"));
138         forAllIter(IOobjectList, vectorFields, iter)
139         {
140             // Read field
141             volVectorField field(*iter(), mesh);
143             // lookup field from dictionary and convert field
144             label unitNumber;
145             if
146             (
147                 foamDataToFluentDict.readIfPresent(field.name(), unitNumber)
148              && unitNumber > 0
149             )
150             {
151                 Info<< "    Converting field " << field.name() << endl;
152                 writeFluentField(field, unitNumber, fluentDataFile);
153             }
154         }
156         Info<< endl;
157     }
159     Info<< "End\n" << endl;
161     return 0;
165 // ************************************************************************* //