FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

createTurbulenceFields.C

Go to the documentation of this file.
00001 /*---------------------------------------------------------------------------*\
00002   =========                 |
00003   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
00004    \\    /   O peration     |
00005     \\  /    A nd           | Copyright (C) 1991-2010 OpenCFD Ltd.
00006      \\/     M anipulation  |
00007 -------------------------------------------------------------------------------
00008 License
00009     This file is part of OpenFOAM.
00010 
00011     OpenFOAM is free software: you can redistribute it and/or modify it
00012     under the terms of the GNU General Public License as published by
00013     the Free Software Foundation, either version 3 of the License, or
00014     (at your option) any later version.
00015 
00016     OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
00017     ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
00018     FITNESS FOR A PARTICULAR PURPOSE.  See the GNU General Public License
00019     for more details.
00020 
00021     You should have received a copy of the GNU General Public License
00022     along with OpenFOAM.  If not, see <http://www.gnu.org/licenses/>.
00023 
00024 Application
00025     createTurbulenceFields
00026 
00027 Description
00028     Creates a full set of turbulence fields.
00029 
00030     - Currently does not output nut and nuTilda
00031 
00032 Source files:
00033     createFields.H
00034 
00035 Usage
00036 
00037     - createTurbulenceFields [OPTIONS]
00038 
00039     @param -noZero \n
00040     Ignore timestep 0.
00041 
00042     @param -constant \n
00043     Include the constant directory.
00044 
00045     @param -time <time>\n
00046     Apply only to specific time.
00047 
00048     @param -latestTime \n
00049     Only apply to latest time step.
00050 
00051     @param -case <dir>\n
00052     Case directory.
00053 
00054     @param -parallel \n
00055     Run in parallel.
00056 
00057     @param -help \n
00058     Display help message.
00059 
00060     @param -doc \n
00061     Display Doxygen API documentation page for this application.
00062 
00063     @param -srcDoc \n
00064     Display Doxygen source documentation page for this application.
00065 
00066 \*---------------------------------------------------------------------------*/
00067 
00068 #include <finiteVolume/fvCFD.H>
00069 #include <incompressibleTransportModels/singlePhaseTransportModel.H>
00070 #include <incompressibleRASModels/RASModel.H>
00071 
00072 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00073 
00074 int main(int argc, char *argv[])
00075 {
00076     timeSelector::addOptions();
00077 
00078     #include <OpenFOAM/setRootCase.H>
00079     #include <OpenFOAM/createTime.H>
00080 
00081     instantList timeDirs = timeSelector::select0(runTime, args);
00082 
00083     #include <OpenFOAM/createMesh.H>
00084     #include "createFields.H"
00085 
00086     forAll(timeDirs, timeI)
00087     {
00088         runTime.setTime(timeDirs[timeI], timeI);
00089 
00090         Info<< "Time = " << runTime.timeName() << endl;
00091 
00092         // Cache the turbulence fields
00093 
00094         Info<< "\nRetrieving field k from turbulence model" << endl;
00095         const volScalarField k = RASModel->k();
00096 
00097         Info<< "\nRetrieving field epsilon from turbulence model" << endl;
00098         const volScalarField epsilon = RASModel->epsilon();
00099 
00100         Info<< "\nRetrieving field R from turbulence model" << endl;
00101         const volSymmTensorField R = RASModel->R();
00102 
00103         // Check availability of tubulence fields
00104 
00105         if (!IOobject("k", runTime.timeName(), mesh).headerOk())
00106         {
00107             Info<< "\nWriting turbulence field k" << endl;
00108             k.write();
00109         }
00110         else
00111         {
00112             Info<< "\nTurbulence k field already exists" << endl;
00113         }
00114 
00115         if (!IOobject("epsilon", runTime.timeName(), mesh).headerOk())
00116         {
00117             Info<< "\nWriting turbulence field epsilon" << endl;
00118             epsilon.write();
00119         }
00120         else
00121         {
00122             Info<< "\nTurbulence epsilon field already exists" << endl;
00123         }
00124 
00125         if (!IOobject("R", runTime.timeName(), mesh).headerOk())
00126         {
00127             Info<< "\nWriting turbulence field R" << endl;
00128             R.write();
00129         }
00130         else
00131         {
00132             Info<< "\nTurbulence R field already exists" << endl;
00133         }
00134 
00135         if (!IOobject("omega", runTime.timeName(), mesh).headerOk())
00136         {
00137             const scalar Cmu = 0.09;
00138 
00139             Info<< "creating omega" << endl;
00140             volScalarField omega
00141             (
00142                 IOobject
00143                 (
00144                     "omega",
00145                     runTime.timeName(),
00146                     mesh
00147                 ),
00148                 epsilon/(Cmu*k),
00149                 epsilon.boundaryField().types()
00150             );
00151             Info<< "\nWriting turbulence field omega" << endl;
00152             omega.write();
00153         }
00154         else
00155         {
00156             Info<< "\nTurbulence omega field already exists" << endl;
00157         }
00158     }
00159 
00160     Info<< "\nEnd\n" << endl;
00161 
00162     return 0;
00163 }
00164 
00165 
00166 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines