1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
25 Translates FOAM data to Fluent format.
27 \*---------------------------------------------------------------------------*/
30 #include "writeFluentFields.H"
32 #include "IOobjectList.H"
34 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
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)
54 runTime.setTime(timeDirs[timeI], timeI);
56 Info<< "Time = " << runTime.timeName() << endl;
58 if (mesh.readUpdate())
60 Info<< " Read new mesh" << endl;
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
72 runTime.caseName() + runTime.timeName() + ".dat"
76 << "(0 \"FOAM to Fluent data File\")" << endl << endl;
78 // Writing number of faces
79 label nFaces = mesh.nFaces();
81 forAll (mesh.boundary(), patchI)
83 nFaces += mesh.boundary()[patchI].size();
87 << "(33 (" << mesh.nCells() << " " << nFaces << " "
88 << mesh.nPoints() << "))" << endl;
90 IOdictionary foamDataToFluentDict
94 "foamDataToFluentDict",
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"));
115 IOobjectList::iterator scalarFieldIter = scalarFields.begin();
116 scalarFieldIter != scalarFields.end();
127 // lookup field from dictionary
128 if (foamDataToFluentDict.found(field.name()))
132 readLabel(foamDataToFluentDict.lookup(field.name()))
138 Info<< " Converting field " << field.name() << endl;
139 writeFluentField(field, unitNumber, fluentDataFile);
145 // Converting volVectorField
146 // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
148 // Search list of objects for volVectorFields
149 IOobjectList vectorFields(objects.lookupClass("volVectorField"));
153 IOobjectList::iterator vectorFieldIter = vectorFields.begin();
154 vectorFieldIter != vectorFields.end();
165 // lookup field from dictionary
166 if (foamDataToFluentDict.found(field.name()))
170 readLabel(foamDataToFluentDict.lookup(field.name()))
176 Info<< " Converting field " << field.name() << endl;
177 writeFluentField(field, unitNumber, fluentDataFile);
185 Info << "End\n" << endl;
191 // ************************************************************************* //