1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | foam-extend: Open Source CFD
4 \\ / O peration | Version: 3.2
5 \\ / A nd | Web: http://www.foam-extend.org
6 \\/ M anipulation | For copyright notice see file Copyright
7 -------------------------------------------------------------------------------
9 This file is part of foam-extend.
11 foam-extend 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 3 of the License, or (at your
14 option) any later version.
16 foam-extend is distributed in the hope that it will be useful, but
17 WITHOUT ANY WARRANTY; without even the implied warranty of
18 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
19 General Public License for more details.
21 You should have received a copy of the GNU General Public License
22 along with foam-extend. If not, see <http://www.gnu.org/licenses/>.
25 CreateTurbulenceFields
28 Creates a full set of turbulence fields.
30 - Currently does not output nut and nuTilda
35 \*---------------------------------------------------------------------------*/
38 #include "incompressible/singlePhaseTransportModel/singlePhaseTransportModel.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)
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())
76 Info<< "\nWriting turbulence field k" << endl;
81 Info<< "\nTurbulence k field already exists" << endl;
84 if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk())
86 Info<< "\nWriting turbulence field epsilon" << endl;
91 Info<< "\nTurbulence epsilon field already exists" << endl;
94 if (!IOobject("R", runTime.timeName(), mesh).headerOk())
96 Info<< "\nWriting turbulence field R" << endl;
101 Info<< "\nTurbulence R field already exists" << endl;
104 if (!IOobject("omega", runTime.timeName(), mesh).headerOk())
106 const scalar Cmu = 0.09;
108 Info<< "creating omega" << endl;
118 epsilon.boundaryField().types()
120 Info<< "\nWriting turbulence field omega" << endl;
125 Info<< "\nTurbulence omega field already exists" << endl;
129 Info<< "\nEnd\n" << endl;
135 // ************************************************************************* //