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

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