initial commit for version 1.6.x patch release
[OpenFOAM-1.6.x.git] / applications / utilities / postProcessing / turbulence / createTurbulenceFields / createTurbulenceFields.C
blob02a8cbd3ab7971c7698e9d75e0379aef45ef2cc9
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2009 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
9     This file is part of OpenFOAM.
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     option) any later version.
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.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 Application
26     CreateTurbulenceFields
28 Description
29     Creates a full set of turbulence fields.
31     - Currently does not output nut and nuTilda
33 Source files:
34     createFields.H
36 \*---------------------------------------------------------------------------*/
38 #include "fvCFD.H"
39 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
40 #include "RASModel.H"
42 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
44 int main(int argc, char *argv[])
46     timeSelector::addOptions();
48     #include "setRootCase.H"
49     #include "createTime.H"
51     instantList timeDirs = timeSelector::select0(runTime, args);
53     #include "createMesh.H"
54     #include "createFields.H"
56     forAll(timeDirs, timeI)
57     {
58         runTime.setTime(timeDirs[timeI], timeI);
60         Info<< "Time = " << runTime.timeName() << endl;
62         // Cache the turbulence fields
64         Info<< "\nRetrieving field k from turbulence model" << endl;
65         const volScalarField k = RASModel->k();
67         Info<< "\nRetrieving field epsilon from turbulence model" << endl;
68         const volScalarField epsilon = RASModel->epsilon();
70         Info<< "\nRetrieving field R from turbulence model" << endl;
71         const volSymmTensorField R = RASModel->R();
73         // Check availability of tubulence fields
75         if (!IOobject("k", runTime.timeName(), mesh).headerOk())
76         {
77             Info<< "\nWriting turbulence field k" << endl;
78             k.write();
79         }
80         else
81         {
82             Info<< "\nTurbulence k field already exists" << endl;
83         }
85         if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk())
86         {
87             Info<< "\nWriting turbulence field epsilon" << endl;
88             epsilon.write();
89         }
90         else
91         {
92             Info<< "\nTurbulence epsilon field already exists" << endl;
93         }
95         if (!IOobject("R", runTime.timeName(), mesh).headerOk())
96         {
97             Info<< "\nWriting turbulence field R" << endl;
98             R.write();
99         }
100         else
101         {
102             Info<< "\nTurbulence R field already exists" << endl;
103         }
105         if (!IOobject("omega", runTime.timeName(), mesh).headerOk())
106         {
107             const scalar Cmu = 0.09;
109             Info<< "creating omega" << endl;
110             volScalarField omega
111             (
112                 IOobject
113                 (
114                     "omega",
115                     runTime.timeName(),
116                     mesh
117                 ),
118                 epsilon/(Cmu*k),
119                 epsilon.boundaryField().types()
120             );
121             Info<< "\nWriting turbulence field omega" << endl;
122             omega.write();
123         }
124         else
125         {
126             Info<< "\nTurbulence omega field already exists" << endl;
127         }
128     }
130     Info<< "\nEnd\n" << endl;
132     return 0;
136 // ************************************************************************* //