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

SpalartAllmaras.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 "SpalartAllmaras.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace incompressible
00034 {
00035 namespace RASModels
00036 {
00037 
00038 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00039 
00040 defineTypeNameAndDebug(SpalartAllmaras, 0);
00041 addToRunTimeSelectionTable(RASModel, SpalartAllmaras, dictionary);
00042 
00043 // * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //
00044 
00045 tmp<volScalarField> SpalartAllmaras::chi() const
00046 {
00047     return nuTilda_/nu();
00048 }
00049 
00050 
00051 tmp<volScalarField> SpalartAllmaras::fv1(const volScalarField& chi) const
00052 {
00053     volScalarField chi3 = pow3(chi);
00054     return chi3/(chi3 + pow3(Cv1_));
00055 }
00056 
00057 
00058 tmp<volScalarField> SpalartAllmaras::fv2
00059 (
00060     const volScalarField& chi,
00061     const volScalarField& fv1
00062 ) const
00063 {
00064     return 1.0/pow3(scalar(1) + chi/Cv2_);
00065 }
00066 
00067 
00068 tmp<volScalarField> SpalartAllmaras::fv3
00069 (
00070     const volScalarField& chi,
00071     const volScalarField& fv1
00072 ) const
00073 {
00074     volScalarField chiByCv2 = (1/Cv2_)*chi;
00075 
00076     return
00077         (scalar(1) + chi*fv1)
00078        *(1/Cv2_)
00079        *(3*(scalar(1) + chiByCv2) + sqr(chiByCv2))
00080        /pow3(scalar(1) + chiByCv2);
00081 }
00082 
00083 
00084 tmp<volScalarField> SpalartAllmaras::fw(const volScalarField& Stilda) const
00085 {
00086     volScalarField r = min
00087     (
00088         nuTilda_
00089        /(
00090            max(Stilda, dimensionedScalar("SMALL", Stilda.dimensions(), SMALL))
00091           *sqr(kappa_*d_)
00092         ),
00093         scalar(10.0)
00094     );
00095     r.boundaryField() == 0.0;
00096 
00097     volScalarField g = r + Cw2_*(pow6(r) - r);
00098 
00099     return g*pow((1.0 + pow6(Cw3_))/(pow6(g) + pow6(Cw3_)), 1.0/6.0);
00100 }
00101 
00102 
00103 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00104 
00105 SpalartAllmaras::SpalartAllmaras
00106 (
00107     const volVectorField& U,
00108     const surfaceScalarField& phi,
00109     transportModel& lamTransportModel
00110 )
00111 :
00112     RASModel(typeName, U, phi, lamTransportModel),
00113 
00114     sigmaNut_
00115     (
00116         dimensioned<scalar>::lookupOrAddToDict
00117         (
00118             "sigmaNut",
00119             coeffDict_,
00120             0.66666
00121         )
00122     ),
00123     kappa_
00124     (
00125         dimensioned<scalar>::lookupOrAddToDict
00126         (
00127             "kappa",
00128             coeffDict_,
00129             0.41
00130         )
00131     ),
00132 
00133     Cb1_
00134     (
00135         dimensioned<scalar>::lookupOrAddToDict
00136         (
00137             "Cb1",
00138             coeffDict_,
00139             0.1355
00140         )
00141     ),
00142     Cb2_
00143     (
00144         dimensioned<scalar>::lookupOrAddToDict
00145         (
00146             "Cb2",
00147             coeffDict_,
00148             0.622
00149         )
00150     ),
00151     Cw1_(Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_),
00152     Cw2_
00153     (
00154         dimensioned<scalar>::lookupOrAddToDict
00155         (
00156             "Cw2",
00157             coeffDict_,
00158             0.3
00159         )
00160     ),
00161     Cw3_
00162     (
00163         dimensioned<scalar>::lookupOrAddToDict
00164         (
00165             "Cw3",
00166             coeffDict_,
00167             2.0
00168         )
00169     ),
00170     Cv1_
00171     (
00172         dimensioned<scalar>::lookupOrAddToDict
00173         (
00174             "Cv1",
00175             coeffDict_,
00176             7.1
00177         )
00178     ),
00179     Cv2_
00180     (
00181         dimensioned<scalar>::lookupOrAddToDict
00182         (
00183             "Cv2",
00184             coeffDict_,
00185             5.0
00186         )
00187     ),
00188 
00189     nuTilda_
00190     (
00191         IOobject
00192         (
00193             "nuTilda",
00194             runTime_.timeName(),
00195             mesh_,
00196             IOobject::MUST_READ,
00197             IOobject::AUTO_WRITE
00198         ),
00199         mesh_
00200     ),
00201 
00202     nut_
00203     (
00204         IOobject
00205         (
00206             "nut",
00207             runTime_.timeName(),
00208             mesh_,
00209             IOobject::MUST_READ,
00210             IOobject::AUTO_WRITE
00211         ),
00212         mesh_
00213     ),
00214 
00215     d_(mesh_)
00216 {
00217     printCoeffs();
00218 }
00219 
00220 
00221 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00222 
00223 tmp<volScalarField> SpalartAllmaras::DnuTildaEff() const
00224 {
00225     return tmp<volScalarField>
00226     (
00227         new volScalarField("DnuTildaEff", nuTilda_/sigmaNut_ + nu())
00228     );
00229 }
00230 
00231 
00232 tmp<volScalarField> SpalartAllmaras::k() const
00233 {
00234     return tmp<volScalarField>
00235     (
00236         new volScalarField
00237         (
00238             IOobject
00239             (
00240                 "k",
00241                 runTime_.timeName(),
00242                 mesh_
00243             ),
00244             mesh_,
00245             dimensionedScalar("0", dimensionSet(0, 2, -2, 0, 0), 0)
00246         )
00247     );
00248 }
00249 
00250 
00251 tmp<volScalarField> SpalartAllmaras::epsilon() const
00252 {
00253     return tmp<volScalarField>
00254     (
00255         new volScalarField
00256         (
00257             IOobject
00258             (
00259                 "epsilon",
00260                 runTime_.timeName(),
00261                 mesh_
00262             ),
00263             mesh_,
00264             dimensionedScalar("0", dimensionSet(0, 2, -3, 0, 0), 0)
00265         )
00266     );
00267 }
00268 
00269 
00270 tmp<volSymmTensorField> SpalartAllmaras::R() const
00271 {
00272     return tmp<volSymmTensorField>
00273     (
00274         new volSymmTensorField
00275         (
00276             IOobject
00277             (
00278                 "R",
00279                 runTime_.timeName(),
00280                 mesh_,
00281                 IOobject::NO_READ,
00282                 IOobject::NO_WRITE
00283             ),
00284             ((2.0/3.0)*I)*k() - nut()*twoSymm(fvc::grad(U_))
00285         )
00286     );
00287 }
00288 
00289 
00290 tmp<volSymmTensorField> SpalartAllmaras::devReff() const
00291 {
00292     return tmp<volSymmTensorField>
00293     (
00294         new volSymmTensorField
00295         (
00296             IOobject
00297             (
00298                 "devRhoReff",
00299                 runTime_.timeName(),
00300                 mesh_,
00301                 IOobject::NO_READ,
00302                 IOobject::NO_WRITE
00303             ),
00304            -nuEff()*dev(twoSymm(fvc::grad(U_)))
00305         )
00306     );
00307 }
00308 
00309 
00310 tmp<fvVectorMatrix> SpalartAllmaras::divDevReff(volVectorField& U) const
00311 {
00312     volScalarField nuEff_ = nuEff();
00313 
00314     return
00315     (
00316       - fvm::laplacian(nuEff_, U)
00317       - fvc::div(nuEff_*dev(fvc::grad(U)().T()))
00318     );
00319 }
00320 
00321 
00322 bool SpalartAllmaras::read()
00323 {
00324     if (RASModel::read())
00325     {
00326         sigmaNut_.readIfPresent(coeffDict());
00327 
00328         Cb1_.readIfPresent(coeffDict());
00329         Cb2_.readIfPresent(coeffDict());
00330         Cw1_ = Cb1_/sqr(kappa_) + (1.0 + Cb2_)/sigmaNut_;
00331         Cw2_.readIfPresent(coeffDict());
00332         Cw3_.readIfPresent(coeffDict());
00333         Cv1_.readIfPresent(coeffDict());
00334         Cv2_.readIfPresent(coeffDict());
00335 
00336         return true;
00337     }
00338     else
00339     {
00340         return false;
00341     }
00342 }
00343 
00344 
00345 void SpalartAllmaras::correct()
00346 {
00347     RASModel::correct();
00348 
00349     if (!turbulence_)
00350     {
00351         return;
00352     }
00353 
00354     if (mesh_.changing())
00355     {
00356         d_.correct();
00357     }
00358 
00359     volScalarField chi = this->chi();
00360     volScalarField fv1 = this->fv1(chi);
00361 
00362     volScalarField Stilda =
00363         fv3(chi, fv1)*::sqrt(2.0)*mag(skew(fvc::grad(U_)))
00364       + fv2(chi, fv1)*nuTilda_/sqr(kappa_*d_);
00365 
00366     tmp<fvScalarMatrix> nuTildaEqn
00367     (
00368         fvm::ddt(nuTilda_)
00369       + fvm::div(phi_, nuTilda_)
00370       - fvm::Sp(fvc::div(phi_), nuTilda_)
00371       - fvm::laplacian(DnuTildaEff(), nuTilda_)
00372       - Cb2_/sigmaNut_*magSqr(fvc::grad(nuTilda_))
00373      ==
00374         Cb1_*Stilda*nuTilda_
00375       - fvm::Sp(Cw1_*fw(Stilda)*nuTilda_/sqr(d_), nuTilda_)
00376     );
00377 
00378     nuTildaEqn().relax();
00379     solve(nuTildaEqn);
00380     bound(nuTilda_, dimensionedScalar("0", nuTilda_.dimensions(), 0.0));
00381     nuTilda_.correctBoundaryConditions();
00382 
00383     nut_.internalField() = fv1*nuTilda_.internalField();
00384     nut_.correctBoundaryConditions();
00385 }
00386 
00387 
00388 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00389 
00390 } // End namespace RASModels
00391 } // End namespace incompressible
00392 } // End namespace Foam
00393 
00394 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines