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

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