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

RASModel.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 \*---------------------------------------------------------------------------*/
00025 
00026 #include "RASModel.H"
00027 #include <finiteVolume/wallFvPatch.H>
00028 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00029 
00030 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00031 
00032 namespace Foam
00033 {
00034 namespace incompressible
00035 {
00036 
00037 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00038 
00039 defineTypeNameAndDebug(RASModel, 0);
00040 defineRunTimeSelectionTable(RASModel, dictionary);
00041 addToRunTimeSelectionTable(turbulenceModel, RASModel, turbulenceModel);
00042 
00043 // * * * * * * * * * * * * * Protected Member Functions  * * * * * * * * * * //
00044 
00045 void RASModel::printCoeffs()
00046 {
00047     if (printCoeffs_)
00048     {
00049         Info<< type() << "Coeffs" << coeffDict_ << endl;
00050     }
00051 }
00052 
00053 
00054 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00055 
00056 RASModel::RASModel
00057 (
00058     const word& type,
00059     const volVectorField& U,
00060     const surfaceScalarField& phi,
00061     transportModel& lamTransportModel
00062 )
00063 :
00064     turbulenceModel(U, phi, lamTransportModel),
00065 
00066     IOdictionary
00067     (
00068         IOobject
00069         (
00070             "RASProperties",
00071             U.time().constant(),
00072             U.db(),
00073             IOobject::MUST_READ,
00074             IOobject::NO_WRITE
00075         )
00076     ),
00077 
00078     turbulence_(lookup("turbulence")),
00079     printCoeffs_(lookupOrDefault<Switch>("printCoeffs", false)),
00080     coeffDict_(subOrEmptyDict(type + "Coeffs")),
00081 
00082     k0_("k0", dimVelocity*dimVelocity, SMALL),
00083     epsilon0_("epsilon0", k0_.dimensions()/dimTime, SMALL),
00084     epsilonSmall_("epsilonSmall", epsilon0_.dimensions(), SMALL),
00085     omega0_("omega0", dimless/dimTime, SMALL),
00086     omegaSmall_("omegaSmall", omega0_.dimensions(), SMALL),
00087 
00088     y_(mesh_)
00089 {
00090     k0_.readIfPresent(*this);
00091     epsilon0_.readIfPresent(*this);
00092     epsilonSmall_.readIfPresent(*this);
00093     omega0_.readIfPresent(*this);
00094     omegaSmall_.readIfPresent(*this);
00095 
00096     // Force the construction of the mesh deltaCoeffs which may be needed
00097     // for the construction of the derived models and BCs
00098     mesh_.deltaCoeffs();
00099 }
00100 
00101 
00102 // * * * * * * * * * * * * * * * * * Selectors * * * * * * * * * * * * * * * //
00103 
00104 autoPtr<RASModel> RASModel::New
00105 (
00106     const volVectorField& U,
00107     const surfaceScalarField& phi,
00108     transportModel& transport
00109 )
00110 {
00111     word modelName;
00112 
00113     // Enclose the creation of the dictionary to ensure it is deleted
00114     // before the turbulenceModel is created otherwise the dictionary is
00115     // entered in the database twice
00116     {
00117         IOdictionary dict
00118         (
00119             IOobject
00120             (
00121                 "RASProperties",
00122                 U.time().constant(),
00123                 U.db(),
00124                 IOobject::MUST_READ,
00125                 IOobject::NO_WRITE
00126             )
00127         );
00128 
00129         dict.lookup("RASModel") >> modelName;
00130     }
00131 
00132     Info<< "Selecting RAS turbulence model " << modelName << endl;
00133 
00134     dictionaryConstructorTable::iterator cstrIter =
00135         dictionaryConstructorTablePtr_->find(modelName);
00136 
00137     if (cstrIter == dictionaryConstructorTablePtr_->end())
00138     {
00139         FatalErrorIn
00140         (
00141             "RASModel::New(const volVectorField&, "
00142             "const surfaceScalarField&, transportModel&)"
00143         )   << "Unknown RASModel type " << modelName
00144             << endl << endl
00145             << "Valid RASModel types are :" << endl
00146             << dictionaryConstructorTablePtr_->sortedToc()
00147             << exit(FatalError);
00148     }
00149 
00150     return autoPtr<RASModel>(cstrIter()(U, phi, transport));
00151 }
00152 
00153 
00154 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00155 
00156 
00157 scalar RASModel::yPlusLam(const scalar kappa, const scalar E) const
00158 {
00159     scalar ypl = 11.0;
00160 
00161     for (int i=0; i<10; i++)
00162     {
00163         ypl = log(max(E*ypl, 1))/kappa;
00164     }
00165 
00166     return ypl;
00167 }
00168 
00169 
00170 tmp<scalarField> RASModel::yPlus(const label patchNo, const scalar Cmu) const
00171 {
00172     const fvPatch& curPatch = mesh_.boundary()[patchNo];
00173 
00174     tmp<scalarField> tYp(new scalarField(curPatch.size()));
00175     scalarField& Yp = tYp();
00176 
00177     if (isA<wallFvPatch>(curPatch))
00178     {
00179         Yp = pow(Cmu, 0.25)
00180             *y_[patchNo]
00181             *sqrt(k()().boundaryField()[patchNo].patchInternalField())
00182             /nu().boundaryField()[patchNo];
00183     }
00184     else
00185     {
00186         WarningIn
00187         (
00188             "tmp<scalarField> RASModel::yPlus(const label patchNo) const"
00189         )   << "Patch " << patchNo << " is not a wall. Returning null field"
00190             << nl << endl;
00191 
00192         Yp.setSize(0);
00193     }
00194 
00195     return tYp;
00196 }
00197 
00198 
00199 void RASModel::correct()
00200 {
00201     turbulenceModel::correct();
00202 
00203     if (turbulence_ && mesh_.changing())
00204     {
00205         y_.correct();
00206     }
00207 }
00208 
00209 
00210 bool RASModel::read()
00211 {
00212     if (regIOobject::read())
00213     {
00214         lookup("turbulence") >> turbulence_;
00215 
00216         if (const dictionary* dictPtr = subDictPtr(type() + "Coeffs"))
00217         {
00218             coeffDict_ <<= *dictPtr;
00219         }
00220 
00221         k0_.readIfPresent(*this);
00222         epsilon0_.readIfPresent(*this);
00223         epsilonSmall_.readIfPresent(*this);
00224         omega0_.readIfPresent(*this);
00225         omegaSmall_.readIfPresent(*this);
00226 
00227         return true;
00228     }
00229     else
00230     {
00231         return false;
00232     }
00233 }
00234 
00235 
00236 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00237 
00238 } // End namespace incompressible
00239 } // End namespace Foam
00240 
00241 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines