FreeFOAM The Cross-Platform CFD Toolkit
createTurbulenceFields.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 | Copyright (C) 1991-2010 OpenCFD Ltd.
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  createTurbulenceFields
26 
27 Description
28  Creates a full set of turbulence fields.
29 
30  - Currently does not output nut and nuTilda
31 
32 Source files:
33  createFields.H
34 
35 Usage
36 
37  - createTurbulenceFields [OPTIONS]
38 
39  @param -noZero \n
40  Ignore timestep 0.
41 
42  @param -constant \n
43  Include the constant directory.
44 
45  @param -time <time>\n
46  Apply only to specific time.
47 
48  @param -latestTime \n
49  Only apply to latest time step.
50 
51  @param -case <dir>\n
52  Case directory.
53 
54  @param -parallel \n
55  Run in parallel.
56 
57  @param -help \n
58  Display help message.
59 
60  @param -doc \n
61  Display Doxygen API documentation page for this application.
62 
63  @param -srcDoc \n
64  Display Doxygen source documentation page for this application.
65 
66 \*---------------------------------------------------------------------------*/
67 
68 #include <finiteVolume/fvCFD.H>
71 
72 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
73 
74 int main(int argc, char *argv[])
75 {
76  timeSelector::addOptions();
77 
78  #include <OpenFOAM/setRootCase.H>
79  #include <OpenFOAM/createTime.H>
80 
81  instantList timeDirs = timeSelector::select0(runTime, args);
82 
83  #include <OpenFOAM/createMesh.H>
84  #include "createFields.H"
85 
86  forAll(timeDirs, timeI)
87  {
88  runTime.setTime(timeDirs[timeI], timeI);
89 
90  Info<< "Time = " << runTime.timeName() << endl;
91 
92  // Cache the turbulence fields
93 
94  Info<< "\nRetrieving field k from turbulence model" << endl;
95  const volScalarField k = RASModel->k();
96 
97  Info<< "\nRetrieving field epsilon from turbulence model" << endl;
98  const volScalarField epsilon = RASModel->epsilon();
99 
100  Info<< "\nRetrieving field R from turbulence model" << endl;
101  const volSymmTensorField R = RASModel->R();
102 
103  // Check availability of tubulence fields
104 
105  if (!IOobject("k", runTime.timeName(), mesh).headerOk())
106  {
107  Info<< "\nWriting turbulence field k" << endl;
108  k.write();
109  }
110  else
111  {
112  Info<< "\nTurbulence k field already exists" << endl;
113  }
114 
115  if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk())
116  {
117  Info<< "\nWriting turbulence field epsilon" << endl;
118  epsilon.write();
119  }
120  else
121  {
122  Info<< "\nTurbulence epsilon field already exists" << endl;
123  }
124 
125  if (!IOobject("R", runTime.timeName(), mesh).headerOk())
126  {
127  Info<< "\nWriting turbulence field R" << endl;
128  R.write();
129  }
130  else
131  {
132  Info<< "\nTurbulence R field already exists" << endl;
133  }
134 
135  if (!IOobject("omega", runTime.timeName(), mesh).headerOk())
136  {
137  const scalar Cmu = 0.09;
138 
139  Info<< "creating omega" << endl;
140  volScalarField omega
141  (
142  IOobject
143  (
144  "omega",
145  runTime.timeName(),
146  mesh
147  ),
148  epsilon/(Cmu*k),
149  epsilon.boundaryField().types()
150  );
151  Info<< "\nWriting turbulence field omega" << endl;
152  omega.write();
153  }
154  else
155  {
156  Info<< "\nTurbulence omega field already exists" << endl;
157  }
158  }
159 
160  Info<< "\nEnd\n" << endl;
161 
162  return 0;
163 }
164 
165 
166 // ************************ vim: set sw=4 sts=4 et: ************************ //