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

LienLeschzinerLowRe.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 "LienLeschzinerLowRe.H"
00027 #include <finiteVolume/wallFvPatch.H>
00028 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00029 
00030 #include <incompressibleRASModels/backwardsCompatibilityWallFunctions.H>
00031 
00032 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00033 
00034 namespace Foam
00035 {
00036 namespace incompressible
00037 {
00038 namespace RASModels
00039 {
00040 
00041 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00042 
00043 defineTypeNameAndDebug(LienLeschzinerLowRe, 0);
00044 addToRunTimeSelectionTable(RASModel, LienLeschzinerLowRe, dictionary);
00045 
00046 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00047 
00048 LienLeschzinerLowRe::LienLeschzinerLowRe
00049 (
00050     const volVectorField& U,
00051     const surfaceScalarField& phi,
00052     transportModel& lamTransportModel
00053 )
00054 :
00055     RASModel(typeName, U, phi, lamTransportModel),
00056 
00057     C1_
00058     (
00059         dimensioned<scalar>::lookupOrAddToDict
00060         (
00061             "C1",
00062             coeffDict_,
00063             1.44
00064         )
00065     ),
00066     C2_
00067     (
00068         dimensioned<scalar>::lookupOrAddToDict
00069         (
00070             "C2",
00071             coeffDict_,
00072             1.92
00073         )
00074     ),
00075     sigmak_
00076     (
00077         dimensioned<scalar>::lookupOrAddToDict
00078         (
00079             "sigmak",
00080             coeffDict_,
00081             1.0
00082         )
00083     ),
00084     sigmaEps_
00085     (
00086         dimensioned<scalar>::lookupOrAddToDict
00087         (
00088             "sigmaEps",
00089             coeffDict_,
00090             1.3
00091         )
00092     ),
00093     Cmu_
00094     (
00095         dimensioned<scalar>::lookupOrAddToDict
00096         (
00097             "Cmu",
00098             coeffDict_,
00099             0.09
00100         )
00101     ),
00102     kappa_
00103     (
00104         dimensioned<scalar>::lookupOrAddToDict
00105         (
00106             "kappa",
00107             coeffDict_,
00108             0.41
00109         )
00110     ),
00111     Am_
00112     (
00113         dimensioned<scalar>::lookupOrAddToDict
00114         (
00115             "Am",
00116             coeffDict_,
00117             0.016
00118         )
00119     ),
00120     Aepsilon_
00121     (
00122         dimensioned<scalar>::lookupOrAddToDict
00123         (
00124             "Aepsilon",
00125             coeffDict_,
00126             0.263
00127         )
00128     ),
00129     Amu_
00130     (
00131         dimensioned<scalar>::lookupOrAddToDict
00132         (
00133             "Amu",
00134             coeffDict_,
00135             0.00222
00136         )
00137     ),
00138 
00139     k_
00140     (
00141         IOobject
00142         (
00143             "k",
00144             runTime_.timeName(),
00145             mesh_,
00146             IOobject::MUST_READ,
00147             IOobject::AUTO_WRITE
00148         ),
00149         mesh_
00150     ),
00151 
00152     epsilon_
00153     (
00154         IOobject
00155         (
00156             "epsilon",
00157             runTime_.timeName(),
00158             mesh_,
00159             IOobject::MUST_READ,
00160             IOobject::AUTO_WRITE
00161         ),
00162         mesh_
00163     ),
00164 
00165     y_(mesh_),
00166 
00167     yStar_(sqrt(k_)*y_/nu() + SMALL),
00168 
00169     nut_
00170     (
00171         IOobject
00172         (
00173             "nut",
00174             runTime_.timeName(),
00175             mesh_,
00176             IOobject::NO_READ,
00177             IOobject::AUTO_WRITE
00178         ),
00179         autoCreateLowReNut("nut", mesh_)
00180     )
00181 {
00182     nut_ = Cmu_*(scalar(1) - exp(-Am_*yStar_))
00183        /(scalar(1) - exp(-Aepsilon_*yStar_) + SMALL)*sqr(k_)
00184        /(epsilon_ + epsilonSmall_);
00185     nut_.correctBoundaryConditions();
00186 
00187     printCoeffs();
00188 }
00189 
00190 
00191 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00192 
00193 tmp<volSymmTensorField> LienLeschzinerLowRe::R() const
00194 {
00195     volTensorField gradU = fvc::grad(U_);
00196 
00197     return tmp<volSymmTensorField>
00198     (
00199         new volSymmTensorField
00200         (
00201             IOobject
00202             (
00203                 "R",
00204                 runTime_.timeName(),
00205                 mesh_,
00206                 IOobject::NO_READ,
00207                 IOobject::NO_WRITE
00208             ),
00209             ((2.0/3.0)*I)*k_ - nut_*twoSymm(gradU),
00210             k_.boundaryField().types()
00211         )
00212     );
00213 }
00214 
00215 
00216 tmp<volSymmTensorField> LienLeschzinerLowRe::devReff() const
00217 {
00218     return tmp<volSymmTensorField>
00219     (
00220         new volSymmTensorField
00221         (
00222             IOobject
00223             (
00224                 "devRhoReff",
00225                 runTime_.timeName(),
00226                 mesh_,
00227                 IOobject::NO_READ,
00228                 IOobject::NO_WRITE
00229             ),
00230            -nuEff()*dev(twoSymm(fvc::grad(U_)))
00231         )
00232     );
00233 }
00234 
00235 
00236 tmp<fvVectorMatrix> LienLeschzinerLowRe::divDevReff(volVectorField& U) const
00237 {
00238     return
00239     (
00240       - fvm::laplacian(nuEff(), U)
00241     //- (fvc::grad(U) & fvc::grad(nuEff()))
00242       - fvc::div(nuEff()*fvc::grad(U)().T())
00243     );
00244 }
00245 
00246 
00247 bool LienLeschzinerLowRe::read()
00248 {
00249     if (RASModel::read())
00250     {
00251         C1_.readIfPresent(coeffDict());
00252         C2_.readIfPresent(coeffDict());
00253         sigmak_.readIfPresent(coeffDict());
00254         sigmaEps_.readIfPresent(coeffDict());
00255         Cmu_.readIfPresent(coeffDict());
00256         kappa_.readIfPresent(coeffDict());
00257         Am_.readIfPresent(coeffDict());
00258         Aepsilon_.readIfPresent(coeffDict());
00259         Amu_.readIfPresent(coeffDict());
00260 
00261         return true;
00262     }
00263     else
00264     {
00265         return false;
00266     }
00267 }
00268 
00269 
00270 void LienLeschzinerLowRe::correct()
00271 {
00272     RASModel::correct();
00273 
00274     if (!turbulence_)
00275     {
00276         return;
00277     }
00278 
00279     if (mesh_.changing())
00280     {
00281         y_.correct();
00282     }
00283 
00284     scalar Cmu75 = pow(Cmu_.value(), 0.75);
00285 
00286     volTensorField gradU = fvc::grad(U_);
00287 
00288     // generation term
00289     volScalarField S2 = symm(gradU) && gradU;
00290 
00291     yStar_ = sqrt(k_)*y_/nu() + SMALL;
00292     volScalarField Rt = sqr(k_)/(nu()*epsilon_);
00293 
00294     volScalarField fMu =
00295         (scalar(1) - exp(-Am_*yStar_))
00296        /(scalar(1) - exp(-Aepsilon_*yStar_) + SMALL);
00297 
00298     volScalarField f2 = scalar(1) - 0.3*exp(-sqr(Rt));
00299 
00300     volScalarField G("RASModel::G", Cmu_*fMu*sqr(k_)/epsilon_*S2);
00301 
00302 
00303     // Dissipation equation
00304     tmp<fvScalarMatrix> epsEqn
00305     (
00306         fvm::ddt(epsilon_)
00307       + fvm::div(phi_, epsilon_)
00308       - fvm::laplacian(DepsilonEff(), epsilon_)
00309       ==
00310         C1_*G*epsilon_/k_
00311         // E-term
00312         + C2_*f2*Cmu75*sqrt(k_)
00313         /(kappa_*y_*(scalar(1) - exp(-Aepsilon_*yStar_)))
00314        *exp(-Amu_*sqr(yStar_))*epsilon_
00315       - fvm::Sp(C2_*f2*epsilon_/k_, epsilon_)
00316     );
00317 
00318     epsEqn().relax();
00319 
00320    #include "LienLeschzinerLowReSetWallDissipation.H"
00321    #include <incompressibleRASModels/wallDissipationI.H>
00322 
00323     solve(epsEqn);
00324     bound(epsilon_, epsilon0_);
00325 
00326 
00327     // Turbulent kinetic energy equation
00328 
00329     tmp<fvScalarMatrix> kEqn
00330     (
00331         fvm::ddt(k_)
00332       + fvm::div(phi_, k_)
00333       - fvm::laplacian(DkEff(), k_)
00334       ==
00335         G
00336       - fvm::Sp(epsilon_/k_, k_)
00337     );
00338 
00339     kEqn().relax();
00340     solve(kEqn);
00341     bound(k_, k0_);
00342 
00343 
00344     // Re-calculate viscosity
00345     nut_ = Cmu_*fMu*sqr(k_)/epsilon_;
00346 }
00347 
00348 
00349 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00350 
00351 } // End namespace RASModels
00352 } // End namespace incompressible
00353 } // End namespace Foam
00354 
00355 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines