foamDataToFluent.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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 Application
27  foamDataToFluent
28 
29 Group
30  grpPostProcessingUtilities
31 
32 Description
33  Translate OpenFOAM data to Fluent format.
34 
35 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
38 #include "writeFluentFields.H"
39 #include "OFstream.H"
40 #include "IOobjectList.H"
41 
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 
44 int main(int argc, char *argv[])
45 {
46  argList::addNote
47  (
48  "Translate OpenFOAM data to Fluent format"
49  );
50  argList::noParallel();
51  timeSelector::addOptions(false); // constant(false), zero(false)
52 
53  #include "setRootCase.H"
54  #include "createTime.H"
55 
56  instantList timeDirs = timeSelector::select0(runTime, args);
57 
58  #include "createNamedMesh.H"
59 
60  // make a directory called proInterface in the case
61  mkDir(runTime.rootPath()/runTime.caseName()/"fluentInterface");
62 
63  forAll(timeDirs, timeI)
64  {
65  runTime.setTime(timeDirs[timeI], timeI);
66 
67  Info<< "Time = " << runTime.timeName() << endl;
68 
69  if (mesh.readUpdate())
70  {
71  Info<< " Read new mesh" << endl;
72  }
73 
74  // make a directory called proInterface in the case
75  mkDir(runTime.rootPath()/runTime.caseName()/"fluentInterface");
76 
77  // open a file for the mesh
78  OFstream fluentDataFile
79  (
80  runTime.rootPath()/
81  runTime.caseName()/
82  "fluentInterface"/
83  runTime.caseName() + runTime.timeName() + ".dat"
84  );
85 
86  fluentDataFile
87  << "(0 \"FOAM to Fluent data File\")" << endl << endl;
88 
89  // Writing number of faces
90  label nFaces = mesh.nFaces();
91 
92  forAll(mesh.boundary(), patchi)
93  {
94  nFaces += mesh.boundary()[patchi].size();
95  }
96 
97  fluentDataFile
98  << "(33 (" << mesh.nCells() << " " << nFaces << " "
99  << mesh.nPoints() << "))" << endl;
100 
101  IOdictionary foamDataToFluentDict
102  (
103  IOobject
104  (
105  "foamDataToFluentDict",
106  runTime.system(),
107  mesh,
108  IOobject::MUST_READ_IF_MODIFIED,
109  IOobject::NO_WRITE
110  )
111  );
112 
113 
114  // Search for list of objects for this time
115  IOobjectList objects(mesh, runTime.timeName());
116 
117 
118  // volScalarField
119  for
120  (
121  const word& fieldName
122  : objects.sortedNames(volScalarField::typeName)
123  )
124  {
125  // Lookup field from dictionary and convert field
126  label unitNumber;
127  if
128  (
129  foamDataToFluentDict.readIfPresent(fieldName, unitNumber)
130  && unitNumber > 0
131  )
132  {
133  // Read field
134  volScalarField field(*(objects[fieldName]), mesh);
135 
136  Info<< " Converting field " << fieldName << nl;
137  writeFluentField(field, unitNumber, fluentDataFile);
138  }
139  }
140 
141 
142  // volVectorField
143  for
144  (
145  const word& fieldName
146  : objects.sortedNames(volVectorField::typeName)
147  )
148  {
149  // Lookup field from dictionary and convert field
150  label unitNumber;
151  if
152  (
153  foamDataToFluentDict.readIfPresent(fieldName, unitNumber)
154  && unitNumber > 0
155  )
156  {
157  // Read field
158  volVectorField field(*(objects[fieldName]), mesh);
159 
160  Info<< " Converting field " << fieldName << nl;
161  writeFluentField(field, unitNumber, fluentDataFile);
162  }
163  }
164 
165  Info<< endl;
166  }
167 
168  Info<< "End\n" << endl;
169 
170  return 0;
171 }
172 
173 
174 // ************************************************************************* //
runTime
engineTime & runTime
Definition: createEngineTime.H:13
IOobjectList.H
Foam::endl
Ostream & endl(Ostream &os)
Definition: Ostream.H:381
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:38
forAll
#define forAll(list, i)
Definition: stdFoam.H:349
OFstream.H
writeFluentFields.H
Foam::Info
messageStream Info
field
rDeltaTY field()
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:53
createNamedMesh.H
Required Variables.
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::volVectorField
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:58
setRootCase.H
Foam::nl
constexpr char nl
Definition: Ostream.H:424
createTime.H
fvCFD.H
args
Foam::argList args(argc, argv)
Foam::mkDir
bool mkDir(const fileName &pathName, mode_t mode=0777)
Definition: POSIX.C:533
Foam::writeFluentField
void writeFluentField(const volScalarField &phi, const label fluentFieldIdentifier, Ostream &stream)