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) 2011-2015 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  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 \*---------------------------------------------------------------------------*/
36 
37 #include "fvCFD.H"
40 
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
42 
43 int main(int argc, char *argv[])
44 {
45  timeSelector::addOptions();
46 
47  argList::addOption
48  (
49  "fields",
50  "wordReList",
51  "specify which turbulence fields (k, epsilon, omega, R) to write"
52  " - eg '(k omega)' or '(R)' or '(.*)'."
53  );
54 
55  #include "setRootCase.H"
56  #include "createTime.H"
57 
58  instantList timeDirs = timeSelector::select0(runTime, args);
59 
60  const bool selectedFields = args.optionFound("fields");
61  wordReList fieldPatterns;
62  if (selectedFields)
63  {
64  fieldPatterns = wordReList(args.optionLookup("fields")());
65  }
66 
67  #include "createMesh.H"
68 
69  forAll(timeDirs, timeI)
70  {
71  runTime.setTime(timeDirs[timeI], timeI);
72 
73  Info<< "Time = " << runTime.timeName() << endl;
74 
75  #include "createFields.H"
76 
77  if (findStrings(fieldPatterns, "k"))
78  {
79  if (!IOobject("k", runTime.timeName(), mesh).headerOk())
80  {
81  Info<< " Writing turbulence field k" << endl;
82  const volScalarField k(RASModel->k());
83  k.write();
84  }
85  else
86  {
87  Info<< " Turbulence k field already exists" << endl;
88  }
89  }
90 
91  if (findStrings(fieldPatterns, "epsilon"))
92  {
93  if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk())
94  {
95  Info<< " Writing turbulence field epsilon" << endl;
96  const volScalarField epsilon(RASModel->epsilon());
97  epsilon.write();
98  }
99  else
100  {
101  Info<< " Turbulence epsilon field already exists" << endl;
102  }
103  }
104 
105  if (findStrings(fieldPatterns, "R"))
106  {
107  if (!IOobject("R", runTime.timeName(), mesh).headerOk())
108  {
109  Info<< " Writing turbulence field R" << endl;
110  const volSymmTensorField R(RASModel->R());
111  R.write();
112  }
113  else
114  {
115  Info<< " Turbulence R field already exists" << endl;
116  }
117  }
118 
119  if (findStrings(fieldPatterns, "omega"))
120  {
121  if (!IOobject("omega", runTime.timeName(), mesh).headerOk())
122  {
123  const scalar Cmu = 0.09;
124 
125  // Assume k and epsilon are available
126  const volScalarField k(RASModel->k());
127  const volScalarField epsilon(RASModel->epsilon());
128 
129  volScalarField omega
130  (
131  IOobject
132  (
133  "omega",
134  runTime.timeName(),
135  mesh
136  ),
137  epsilon/(Cmu*k),
138  epsilon.boundaryField().types()
139  );
140 
141  Info<< " Writing turbulence field omega" << endl;
142  omega.write();
143  }
144  else
145  {
146  Info<< " Turbulence omega field already exists" << endl;
147  }
148  }
149  }
150 
151  Info<< "\nEnd\n" << endl;
152 
153  return 0;
154 }
155 
156 
157 // ************************************************************************* //
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
singlePhaseTransportModel.H
turbulentTransportModel.H
Foam::wordReList
List< wordRe > wordReList
A List of wordRe (word or regular expression)
Definition: wordReList.H:50
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::instantList
List< instant > instantList
List of instants.
Definition: instantList.H:42
R
#define R(A, B, C, D, E, F, K, M)
Foam::volSymmTensorField
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:58
Foam::Ostream::write
virtual Ostream & write(const token &)=0
Write next token to stream.
Foam::Info
messageStream Info
main
int main(int argc, char *argv[])
Definition: postCalc.C:54
Foam::volScalarField
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
setRootCase.H
timeDirs
static instantList timeDirs
Definition: globalFoam.H:44
Foam::findStrings
bool findStrings(const wordReListMatcher &matcher, const std::string &str)
Return true if string matches one of the regular expressions.
Definition: stringListOps.H:52
epsilon
epsilon
Definition: createTDFields.H:56
RASModel
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);singlePhaseTransportModel laminarTransport(U, phi);autoPtr< incompressible::RASModel > RASModel(incompressible::New< incompressible::RASModel >(U, phi, laminarTransport))
createMesh.H
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
createTime.H
Foam::argList::optionFound
bool optionFound(const word &opt) const
Return true if the named option is found.
Definition: argListI.H:108
fvCFD.H
args
Foam::argList args(argc, argv)
Foam::argList::optionLookup
IStringStream optionLookup(const word &opt) const
Return an IStringStream from the named option.
Definition: argListI.H:114