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

execFlowFunctionObjects.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 Global
00025     execFlowFunctionObjects
00026 
00027 Application
00028     execFlowFunctionObjects
00029 
00030 Description
00031     Execute the set of functionObjects specified in the selected dictionary
00032     (which defaults to system/controlDict) for the selected set of times.
00033 
00034     The flow (p-U) and optionally turbulence fields are available for the
00035     function objects to operate on allowing forces and other related properties
00036     to be calculated in addition to cutting planes etc.
00037 
00038 Usage
00039 
00040     - execFlowFunctionObjects [OPTIONS]
00041 
00042     @param -noWrite \n
00043     Suppress output to files.
00044 
00045     @param -dict <dictionary name>\n
00046     Use named dictionary instead of system/controlDict.
00047 
00048     @param -noZero \n
00049     Ignore timestep 0.
00050 
00051     @param -constant \n
00052     Include the constant directory.
00053 
00054     @param -time <time>\n
00055     Apply only to specific time.
00056 
00057     @param -latestTime \n
00058     Only apply to latest time step.
00059 
00060     @param -case <dir>\n
00061     Case directory.
00062 
00063     @param -parallel \n
00064     Run in parallel.
00065 
00066     @param -help \n
00067     Display help message.
00068 
00069     @param -doc \n
00070     Display Doxygen API documentation page for this application.
00071 
00072     @param -srcDoc \n
00073     Display Doxygen source documentation page for this application.
00074 
00075 \*---------------------------------------------------------------------------*/
00076 
00077 #include <postCalc/calc.H>
00078 
00079 #include <incompressibleTransportModels/singlePhaseTransportModel.H>
00080 
00081 #include <incompressibleRASModels/RASModel.H>
00082 #include <incompressibleLESModels/LESModel.H>
00083 
00084 #include <basicThermophysicalModels/basicPsiThermo.H>
00085 #include <compressibleRASModels/RASModel.H>
00086 #include <compressibleLESModels/LESModel.H>
00087 
00088 
00089 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00090 
00091 namespace Foam
00092 {
00093     void execFlowFunctionObjects(const argList& args, const Time& runTime)
00094     {
00095         if (args.optionFound("dict"))
00096         {
00097             IOdictionary dict
00098             (
00099                 IOobject
00100                 (
00101                     args.option("dict"),
00102                     runTime.system(),
00103                     runTime,
00104                     IOobject::MUST_READ
00105                 )
00106             );
00107 
00108             functionObjectList fol(runTime, dict);
00109             fol.start();
00110             fol.execute();
00111         }
00112         else
00113         {
00114             functionObjectList fol(runTime);
00115             fol.start();
00116             fol.execute();
00117         }
00118     }
00119 }
00120 
00121 
00122 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00123 
00124 void Foam::calc(const argList& args, const Time& runTime, const fvMesh& mesh)
00125 {
00126     Info<< "    Reading phi" << endl;
00127     surfaceScalarField phi
00128     (
00129         IOobject
00130         (
00131             "phi",
00132             runTime.timeName(),
00133             mesh,
00134             IOobject::MUST_READ
00135         ),
00136         mesh
00137     );
00138 
00139     Info<< "    Reading U" << endl;
00140     volVectorField U
00141     (
00142         IOobject
00143         (
00144             "U",
00145             runTime.timeName(),
00146             mesh,
00147             IOobject::MUST_READ
00148         ),
00149         mesh
00150     );
00151 
00152     Info<< "    Reading p" << endl;
00153     volScalarField p
00154     (
00155         IOobject
00156         (
00157             "p",
00158             runTime.timeName(),
00159             mesh,
00160             IOobject::MUST_READ
00161         ),
00162         mesh
00163     );
00164 
00165     if (phi.dimensions() == dimensionSet(0, 3, -1, 0, 0))
00166     {
00167         IOobject RASPropertiesHeader
00168         (
00169             "RASProperties",
00170             runTime.constant(),
00171             mesh,
00172             IOobject::MUST_READ,
00173             IOobject::NO_WRITE,
00174             false
00175         );
00176 
00177         IOobject LESPropertiesHeader
00178         (
00179             "LESProperties",
00180             runTime.constant(),
00181             mesh,
00182             IOobject::MUST_READ,
00183             IOobject::NO_WRITE,
00184             false
00185         );
00186 
00187         singlePhaseTransportModel laminarTransport(U, phi);
00188 
00189         if (RASPropertiesHeader.headerOk())
00190         {
00191             IOdictionary RASProperties(RASPropertiesHeader);
00192 
00193             autoPtr<incompressible::RASModel> RASModel
00194             (
00195                 incompressible::RASModel::New
00196                 (
00197                     U,
00198                     phi,
00199                     laminarTransport
00200                 )
00201             );
00202             execFlowFunctionObjects(args, runTime);
00203         }
00204         else if (LESPropertiesHeader.headerOk())
00205         {
00206             IOdictionary LESProperties(LESPropertiesHeader);
00207 
00208             autoPtr<incompressible::LESModel> sgsModel
00209             (
00210                 incompressible::LESModel::New(U, phi, laminarTransport)
00211             );
00212 
00213             execFlowFunctionObjects(args, runTime);
00214         }
00215         else
00216         {
00217             IOdictionary transportProperties
00218             (
00219                 IOobject
00220                 (
00221                     "transportProperties",
00222                     runTime.constant(),
00223                     mesh,
00224                     IOobject::MUST_READ,
00225                     IOobject::NO_WRITE
00226                 )
00227             );
00228 
00229             dimensionedScalar nu(transportProperties.lookup("nu"));
00230 
00231             execFlowFunctionObjects(args, runTime);
00232         }
00233     }
00234     else if (phi.dimensions() == dimensionSet(1, 0, -1, 0, 0))
00235     {
00236         autoPtr<basicPsiThermo> thermo(basicPsiThermo::New(mesh));
00237 
00238         volScalarField rho
00239         (
00240             IOobject
00241             (
00242                 "rho",
00243                 runTime.timeName(),
00244                 mesh
00245             ),
00246             thermo->rho()
00247         );
00248 
00249         IOobject RASPropertiesHeader
00250         (
00251             "RASProperties",
00252             runTime.constant(),
00253             mesh,
00254             IOobject::MUST_READ,
00255             IOobject::NO_WRITE,
00256             false
00257         );
00258 
00259         IOobject LESPropertiesHeader
00260         (
00261             "LESProperties",
00262             runTime.constant(),
00263             mesh,
00264             IOobject::MUST_READ,
00265             IOobject::NO_WRITE,
00266             false
00267         );
00268 
00269         if (RASPropertiesHeader.headerOk())
00270         {
00271             IOdictionary RASProperties(RASPropertiesHeader);
00272 
00273             autoPtr<compressible::RASModel> RASModel
00274             (
00275                 compressible::RASModel::New
00276                 (
00277                     rho,
00278                     U,
00279                     phi,
00280                     thermo()
00281                 )
00282             );
00283 
00284             execFlowFunctionObjects(args, runTime);
00285         }
00286         else if (LESPropertiesHeader.headerOk())
00287         {
00288             IOdictionary LESProperties(LESPropertiesHeader);
00289 
00290             autoPtr<compressible::LESModel> sgsModel
00291             (
00292                 compressible::LESModel::New(rho, U, phi, thermo())
00293             );
00294 
00295             execFlowFunctionObjects(args, runTime);
00296         }
00297         else
00298         {
00299             IOdictionary transportProperties
00300             (
00301                 IOobject
00302                 (
00303                     "transportProperties",
00304                     runTime.constant(),
00305                     mesh,
00306                     IOobject::MUST_READ,
00307                     IOobject::NO_WRITE
00308                 )
00309             );
00310 
00311             dimensionedScalar mu(transportProperties.lookup("mu"));
00312 
00313             execFlowFunctionObjects(args, runTime);
00314         }
00315     }
00316     else
00317     {
00318         FatalErrorIn(args.executable())
00319             << "Incorrect dimensions of phi: " << phi.dimensions()
00320             << nl << exit(FatalError);
00321     }
00322 }
00323 
00324 
00325 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines