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

realizableKE.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 "realizableKE.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 #include <compressibleRASModels/backwardsCompatibilityWallFunctions.H>
00030 
00031 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035 namespace compressible
00036 {
00037 namespace RASModels
00038 {
00039 
00040 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
00041 
00042 defineTypeNameAndDebug(realizableKE, 0);
00043 addToRunTimeSelectionTable(RASModel, realizableKE, dictionary);
00044 
00045 // * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //
00046 
00047 tmp<volScalarField> realizableKE::rCmu
00048 (
00049     const volTensorField& gradU,
00050     const volScalarField& S2,
00051     const volScalarField& magS
00052 )
00053 {
00054     tmp<volSymmTensorField> tS = dev(symm(gradU));
00055     const volSymmTensorField& S = tS();
00056 
00057     volScalarField W =
00058         (2*sqrt(2.0))*((S&S)&&S)
00059        /(
00060             magS*S2
00061           + dimensionedScalar("small", dimensionSet(0, 0, -3, 0, 0), SMALL)
00062         );
00063 
00064     tS.clear();
00065 
00066     volScalarField phis =
00067         (1.0/3.0)*acos(min(max(sqrt(6.0)*W, -scalar(1)), scalar(1)));
00068     volScalarField As = sqrt(6.0)*cos(phis);
00069     volScalarField Us = sqrt(S2/2.0 + magSqr(skew(gradU)));
00070 
00071     return 1.0/(A0_ + As*Us*k_/(epsilon_ + epsilonSmall_));
00072 }
00073 
00074 
00075 tmp<volScalarField> realizableKE::rCmu
00076 (
00077     const volTensorField& gradU
00078 )
00079 {
00080     volScalarField S2 = 2*magSqr(dev(symm(gradU)));
00081     volScalarField magS = sqrt(S2);
00082     return rCmu(gradU, S2, magS);
00083 }
00084 
00085 
00086 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00087 
00088 realizableKE::realizableKE
00089 (
00090     const volScalarField& rho,
00091     const volVectorField& U,
00092     const surfaceScalarField& phi,
00093     const basicThermo& thermophysicalModel
00094 )
00095 :
00096     RASModel(typeName, rho, U, phi, thermophysicalModel),
00097 
00098     Cmu_
00099     (
00100         dimensioned<scalar>::lookupOrAddToDict
00101         (
00102             "Cmu",
00103             coeffDict_,
00104             0.09
00105         )
00106     ),
00107     A0_
00108     (
00109         dimensioned<scalar>::lookupOrAddToDict
00110         (
00111             "A0",
00112             coeffDict_,
00113             4.0
00114         )
00115     ),
00116     C2_
00117     (
00118         dimensioned<scalar>::lookupOrAddToDict
00119         (
00120             "C2",
00121             coeffDict_,
00122             1.9
00123         )
00124     ),
00125     sigmak_
00126     (
00127         dimensioned<scalar>::lookupOrAddToDict
00128         (
00129             "sigmak",
00130             coeffDict_,
00131             1.0
00132         )
00133     ),
00134     sigmaEps_
00135     (
00136         dimensioned<scalar>::lookupOrAddToDict
00137         (
00138             "sigmaEps",
00139             coeffDict_,
00140             1.2
00141         )
00142     ),
00143     Prt_
00144     (
00145         dimensioned<scalar>::lookupOrAddToDict
00146         (
00147             "Prt",
00148             coeffDict_,
00149             1.0
00150         )
00151     ),
00152 
00153     k_
00154     (
00155         IOobject
00156         (
00157             "k",
00158             runTime_.timeName(),
00159             mesh_,
00160             IOobject::NO_READ,
00161             IOobject::AUTO_WRITE
00162         ),
00163         autoCreateK("k", mesh_)
00164     ),
00165     epsilon_
00166     (
00167         IOobject
00168         (
00169             "epsilon",
00170             runTime_.timeName(),
00171             mesh_,
00172             IOobject::NO_READ,
00173             IOobject::AUTO_WRITE
00174         ),
00175         autoCreateEpsilon("epsilon", mesh_)
00176     ),
00177     mut_
00178     (
00179         IOobject
00180         (
00181             "mut",
00182             runTime_.timeName(),
00183             mesh_,
00184             IOobject::NO_READ,
00185             IOobject::AUTO_WRITE
00186         ),
00187         autoCreateMut("mut", mesh_)
00188     ),
00189     alphat_
00190     (
00191         IOobject
00192         (
00193             "alphat",
00194             runTime_.timeName(),
00195             mesh_,
00196             IOobject::NO_READ,
00197             IOobject::AUTO_WRITE
00198         ),
00199         autoCreateAlphat("alphat", mesh_)
00200     )
00201 {
00202     bound(k_, k0_);
00203     bound(epsilon_, epsilon0_);
00204 
00205     mut_ = rCmu(fvc::grad(U_))*rho_*sqr(k_)/(epsilon_ + epsilonSmall_);
00206     mut_.correctBoundaryConditions();
00207 
00208     alphat_ = mut_/Prt_;
00209     alphat_.correctBoundaryConditions();
00210 
00211     printCoeffs();
00212 }
00213 
00214 
00215 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00216 
00217 tmp<volSymmTensorField> realizableKE::R() const
00218 {
00219     return tmp<volSymmTensorField>
00220     (
00221         new volSymmTensorField
00222         (
00223             IOobject
00224             (
00225                 "R",
00226                 runTime_.timeName(),
00227                 mesh_,
00228                 IOobject::NO_READ,
00229                 IOobject::NO_WRITE
00230             ),
00231             ((2.0/3.0)*I)*k_ - (mut_/rho_)*dev(twoSymm(fvc::grad(U_))),
00232             k_.boundaryField().types()
00233         )
00234     );
00235 }
00236 
00237 
00238 tmp<volSymmTensorField> realizableKE::devRhoReff() const
00239 {
00240     return tmp<volSymmTensorField>
00241     (
00242         new volSymmTensorField
00243         (
00244             IOobject
00245             (
00246                 "devRhoReff",
00247                 runTime_.timeName(),
00248                 mesh_,
00249                 IOobject::NO_READ,
00250                 IOobject::NO_WRITE
00251             ),
00252            -muEff()*dev(twoSymm(fvc::grad(U_)))
00253         )
00254     );
00255 }
00256 
00257 
00258 tmp<fvVectorMatrix> realizableKE::divDevRhoReff(volVectorField& U) const
00259 {
00260     return
00261     (
00262       - fvm::laplacian(muEff(), U) - fvc::div(muEff()*dev2(fvc::grad(U)().T()))
00263     );
00264 }
00265 
00266 
00267 bool realizableKE::read()
00268 {
00269     if (RASModel::read())
00270     {
00271         Cmu_.readIfPresent(coeffDict());
00272         A0_.readIfPresent(coeffDict());
00273         C2_.readIfPresent(coeffDict());
00274         sigmak_.readIfPresent(coeffDict());
00275         sigmaEps_.readIfPresent(coeffDict());
00276         Prt_.readIfPresent(coeffDict());
00277 
00278         return true;
00279     }
00280     else
00281     {
00282         return false;
00283     }
00284 }
00285 
00286 
00287 void realizableKE::correct()
00288 {
00289     if (!turbulence_)
00290     {
00291         // Re-calculate viscosity
00292         mut_ = rCmu(fvc::grad(U_))*rho_*sqr(k_)/epsilon_;
00293         mut_.correctBoundaryConditions();
00294 
00295         // Re-calculate thermal diffusivity
00296         alphat_ = mut_/Prt_;
00297         alphat_.correctBoundaryConditions();
00298 
00299         return;
00300     }
00301 
00302     RASModel::correct();
00303 
00304     volScalarField divU = fvc::div(phi_/fvc::interpolate(rho_));
00305 
00306     if (mesh_.moving())
00307     {
00308         divU += fvc::div(mesh_.phi());
00309     }
00310 
00311     volTensorField gradU = fvc::grad(U_);
00312     volScalarField S2 = 2*magSqr(dev(symm(gradU)));
00313     volScalarField magS = sqrt(S2);
00314 
00315     volScalarField eta = magS*k_/epsilon_;
00316     volScalarField C1 = max(eta/(scalar(5) + eta), scalar(0.43));
00317 
00318     volScalarField G("RASModel::G", mut_*(gradU && dev(twoSymm(gradU))));
00319 
00320     // Update espsilon and G at the wall
00321     epsilon_.boundaryField().updateCoeffs();
00322 
00323     // Dissipation equation
00324     tmp<fvScalarMatrix> epsEqn
00325     (
00326         fvm::ddt(rho_, epsilon_)
00327       + fvm::div(phi_, epsilon_)
00328       - fvm::laplacian(DepsilonEff(), epsilon_)
00329      ==
00330         C1*rho_*magS*epsilon_
00331       - fvm::Sp
00332         (
00333             C2_*rho_*epsilon_/(k_ + sqrt((mu()/rho_)*epsilon_)),
00334             epsilon_
00335         )
00336     );
00337 
00338     epsEqn().relax();
00339 
00340     epsEqn().boundaryManipulate(epsilon_.boundaryField());
00341 
00342     solve(epsEqn);
00343     bound(epsilon_, epsilon0_);
00344 
00345 
00346     // Turbulent kinetic energy equation
00347 
00348     tmp<fvScalarMatrix> kEqn
00349     (
00350         fvm::ddt(rho_, k_)
00351       + fvm::div(phi_, k_)
00352       - fvm::laplacian(DkEff(), k_)
00353      ==
00354         G - fvm::SuSp(2.0/3.0*rho_*divU, k_)
00355       - fvm::Sp(rho_*(epsilon_)/k_, k_)
00356     );
00357 
00358     kEqn().relax();
00359     solve(kEqn);
00360     bound(k_, k0_);
00361 
00362     // Re-calculate viscosity
00363     mut_ = rCmu(gradU, S2, magS)*rho_*sqr(k_)/epsilon_;
00364     mut_.correctBoundaryConditions();
00365 
00366     // Re-calculate thermal diffusivity
00367     alphat_ = mut_/Prt_;
00368     alphat_.correctBoundaryConditions();
00369 }
00370 
00371 
00372 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00373 
00374 } // End namespace RASModels
00375 } // End namespace compressible
00376 } // End namespace Foam
00377 
00378 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines