ENH: autoLayerDriver: better layering information message
[OpenFOAM-2.0.x.git] / applications / utilities / postProcessing / turbulence / createTurbulenceFields / createTurbulenceFields.C
blob928f0abf5e0e452251ab4a514d2fd104d400eb4a
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 2011 OpenFOAM Foundation
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
13     the Free Software Foundation, either version 3 of the License, or
14     (at your 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, see <http://www.gnu.org/licenses/>.
24 Application
25     CreateTurbulenceFields
27 Description
28     Creates a full set of turbulence fields.
30     - Currently does not output nut and nuTilda
32 Source files:
33     createFields.H
35 \*---------------------------------------------------------------------------*/
37 #include "fvCFD.H"
38 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.H"
39 #include "RASModel.H"
41 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
43 int main(int argc, char *argv[])
45     timeSelector::addOptions();
47     #include "setRootCase.H"
48     #include "createTime.H"
50     instantList timeDirs = timeSelector::select0(runTime, args);
52     #include "createMesh.H"
53     #include "createFields.H"
55     forAll(timeDirs, timeI)
56     {
57         runTime.setTime(timeDirs[timeI], timeI);
59         Info<< "Time = " << runTime.timeName() << endl;
61         // Cache the turbulence fields
63         Info<< "\nRetrieving field k from turbulence model" << endl;
64         const volScalarField k(RASModel->k());
66         Info<< "\nRetrieving field epsilon from turbulence model" << endl;
67         const volScalarField epsilon(RASModel->epsilon());
69         Info<< "\nRetrieving field R from turbulence model" << endl;
70         const volSymmTensorField R(RASModel->R());
72         // Check availability of tubulence fields
74         if (!IOobject("k", runTime.timeName(), mesh).headerOk())
75         {
76             Info<< "\nWriting turbulence field k" << endl;
77             k.write();
78         }
79         else
80         {
81             Info<< "\nTurbulence k field already exists" << endl;
82         }
84         if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk())
85         {
86             Info<< "\nWriting turbulence field epsilon" << endl;
87             epsilon.write();
88         }
89         else
90         {
91             Info<< "\nTurbulence epsilon field already exists" << endl;
92         }
94         if (!IOobject("R", runTime.timeName(), mesh).headerOk())
95         {
96             Info<< "\nWriting turbulence field R" << endl;
97             R.write();
98         }
99         else
100         {
101             Info<< "\nTurbulence R field already exists" << endl;
102         }
104         if (!IOobject("omega", runTime.timeName(), mesh).headerOk())
105         {
106             const scalar Cmu = 0.09;
108             Info<< "creating omega" << endl;
109             volScalarField omega
110             (
111                 IOobject
112                 (
113                     "omega",
114                     runTime.timeName(),
115                     mesh
116                 ),
117                 epsilon/(Cmu*k),
118                 epsilon.boundaryField().types()
119             );
120             Info<< "\nWriting turbulence field omega" << endl;
121             omega.write();
122         }
123         else
124         {
125             Info<< "\nTurbulence omega field already exists" << endl;
126         }
127     }
129     Info<< "\nEnd\n" << endl;
131     return 0;
135 // ************************************************************************* //