foamDataToFluent.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 Application
25  foamDataToFluent
26 
27 Description
28  Translates OpenFOAM data to Fluent format.
29 
30 \*---------------------------------------------------------------------------*/
31 
32 #include "argList.H"
33 #include "timeSelector.H"
34 #include "writeFluentFields.H"
35 #include "OFstream.H"
36 #include "IOobjectList.H"
37 
38 using namespace Foam;
39 
40 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 
42 int main(int argc, char *argv[])
43 {
45  timeSelector::addOptions(false); // no constant
46 
47  #include "setRootCase.H"
48  #include "createTime.H"
49 
51 
52  #include "createMeshNoChangers.H"
53 
54  // make a directory called proInterface in the case
55  mkDir(runTime.rootPath()/runTime.caseName()/"fluentInterface");
56 
57  forAll(timeDirs, timeI)
58  {
59  runTime.setTime(timeDirs[timeI], timeI);
60 
61  Info<< "Time = " << runTime.userTimeName() << endl;
62 
63  if (mesh.readUpdate())
64  {
65  Info<< " Read new mesh" << endl;
66  }
67 
68  // make a directory called proInterface in the case
69  mkDir(runTime.rootPath()/runTime.caseName()/"fluentInterface");
70 
71  // open a file for the mesh
72  OFstream fluentDataFile
73  (
74  runTime.rootPath()/
75  runTime.caseName()/
76  "fluentInterface"/
77  runTime.caseName() + runTime.name() + ".dat"
78  );
79 
80  fluentDataFile
81  << "(0 \"FOAM to Fluent data File\")" << endl << endl;
82 
83  // Writing number of faces
84  label nFaces = mesh.nFaces();
85 
86  forAll(mesh.boundary(), patchi)
87  {
88  nFaces += mesh.boundary()[patchi].size();
89  }
90 
91  fluentDataFile
92  << "(33 (" << mesh.nCells() << " " << nFaces << " "
93  << mesh.nPoints() << "))" << endl;
94 
95  IOdictionary foamDataToFluentDict
96  (
97  IOobject
98  (
99  "foamDataToFluentDict",
100  runTime.system(),
101  mesh,
104  )
105  );
106 
107 
108  // Search for list of objects for this time
109  IOobjectList objects(mesh, runTime.name());
110 
111 
112  // Converting volScalarField
113  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
114 
115  // Search list of objects for volScalarFields
116  IOobjectList scalarFields(objects.lookupClass("volScalarField"));
117 
118  forAllIter(IOobjectList, scalarFields, iter)
119  {
120  // Read field
121  volScalarField field(*iter(), mesh);
122 
123  // lookup field from dictionary and convert field
124  label unitNumber;
125  if
126  (
127  foamDataToFluentDict.readIfPresent(field.name(), unitNumber)
128  && unitNumber > 0
129  )
130  {
131  Info<< " Converting field " << field.name() << endl;
132  writeFluentField(field, unitNumber, fluentDataFile);
133  }
134  }
135 
136 
137  // Converting volVectorField
138  // ~~~~~~~~~~~~~~~~~~~~~~~~~~~~
139 
140  // Search list of objects for volVectorFields
141  IOobjectList vectorFields(objects.lookupClass("volVectorField"));
142 
143  forAllIter(IOobjectList, vectorFields, iter)
144  {
145  // Read field
146  volVectorField field(*iter(), mesh);
147 
148  // lookup field from dictionary and convert field
149  label unitNumber;
150  if
151  (
152  foamDataToFluentDict.readIfPresent(field.name(), unitNumber)
153  && unitNumber > 0
154  )
155  {
156  Info<< " Converting field " << field.name() << endl;
157  writeFluentField(field, unitNumber, fluentDataFile);
158  }
159  }
160 
161  Info<< endl;
162  }
163 
164  Info<< "End\n" << endl;
165 
166  return 0;
167 }
168 
169 
170 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
Generic GeometricField class.
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
List of IOobjects with searching and retrieving facilities.
Definition: IOobjectList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
@ MUST_READ_IF_MODIFIED
Definition: IOobject.H:119
Output to file stream.
Definition: OFstream.H:86
static void noParallel()
Remove the parallel options.
Definition: argList.C:175
static void addOptions(const bool constant=true, const bool withZero=false)
Add the options handled by timeSelector to argList::validOptions.
Definition: timeSelector.C:114
static instantList select0(Time &runTime, const argList &args)
Return the set of times selected based on the argList options.
Definition: timeSelector.C:252
int main(int argc, char *argv[])
Definition: financialFoam.C:44
label patchi
static instantList timeDirs
Definition: globalFoam.H:44
Namespace for OpenFOAM.
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
bool mkDir(const fileName &, mode_t=0777)
Make a directory and return an error if it could not be created.
Definition: POSIX.C:290
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
messageStream Info
void writeFluentField(const volScalarField &phi, const label fluentFieldIdentifier, Ostream &stream)
objects
Foam::argList args(argc, argv)