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 <incompressibleRASModels/backwardsCompatibilityWallFunctions.H>
00030 
00031 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00032 
00033 namespace Foam
00034 {
00035 namespace incompressible
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 volVectorField& U,
00050     const surfaceScalarField& phi,
00051     transportModel& lamTransportModel
00052 )
00053 :
00054     RASModel(typeName, U, phi, lamTransportModel),
00055 
00056     Cmu_
00057     (
00058         dimensioned<scalar>::lookupOrAddToDict
00059         (
00060             "Cmu",
00061             coeffDict_,
00062             0.09
00063         )
00064     ),
00065     C1_
00066     (
00067         dimensioned<scalar>::lookupOrAddToDict
00068         (
00069             "C1",
00070             coeffDict_,
00071             1.44
00072         )
00073     ),
00074     C2_
00075     (
00076         dimensioned<scalar>::lookupOrAddToDict
00077         (
00078             "C2",
00079             coeffDict_,
00080             1.92
00081         )
00082     ),
00083     sigmaEps_
00084     (
00085         dimensioned<scalar>::lookupOrAddToDict
00086         (
00087             "sigmaEps",
00088             coeffDict_,
00089             1.3
00090         )
00091     ),
00092 
00093     k_
00094     (
00095         IOobject
00096         (
00097             "k",
00098             runTime_.timeName(),
00099             mesh_,
00100             IOobject::NO_READ,
00101             IOobject::AUTO_WRITE
00102         ),
00103         autoCreateK("k", mesh_)
00104     ),
00105     epsilon_
00106     (
00107         IOobject
00108         (
00109             "epsilon",
00110             runTime_.timeName(),
00111             mesh_,
00112             IOobject::NO_READ,
00113             IOobject::AUTO_WRITE
00114         ),
00115         autoCreateEpsilon("epsilon", mesh_)
00116     ),
00117     nut_
00118     (
00119         IOobject
00120         (
00121             "nut",
00122             runTime_.timeName(),
00123             mesh_,
00124             IOobject::NO_READ,
00125             IOobject::AUTO_WRITE
00126         ),
00127         autoCreateNut("nut", mesh_)
00128     )
00129 {
00130     nut_ = Cmu_*sqr(k_)/(epsilon_ + epsilonSmall_);
00131     nut_.correctBoundaryConditions();
00132 
00133     printCoeffs();
00134 }
00135 
00136 
00137 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00138 
00139 tmp<volSymmTensorField> kEpsilon::R() const
00140 {
00141     return tmp<volSymmTensorField>
00142     (
00143         new volSymmTensorField
00144         (
00145             IOobject
00146             (
00147                 "R",
00148                 runTime_.timeName(),
00149                 mesh_,
00150                 IOobject::NO_READ,
00151                 IOobject::NO_WRITE
00152             ),
00153             ((2.0/3.0)*I)*k_ - nut_*twoSymm(fvc::grad(U_)),
00154             k_.boundaryField().types()
00155         )
00156     );
00157 }
00158 
00159 
00160 tmp<volSymmTensorField> kEpsilon::devReff() const
00161 {
00162     return tmp<volSymmTensorField>
00163     (
00164         new volSymmTensorField
00165         (
00166             IOobject
00167             (
00168                 "devRhoReff",
00169                 runTime_.timeName(),
00170                 mesh_,
00171                 IOobject::NO_READ,
00172                 IOobject::NO_WRITE
00173             ),
00174            -nuEff()*dev(twoSymm(fvc::grad(U_)))
00175         )
00176     );
00177 }
00178 
00179 
00180 tmp<fvVectorMatrix> kEpsilon::divDevReff(volVectorField& U) const
00181 {
00182     return
00183     (
00184       - fvm::laplacian(nuEff(), U)
00185       - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00186     );
00187 }
00188 
00189 
00190 bool kEpsilon::read()
00191 {
00192     if (RASModel::read())
00193     {
00194         Cmu_.readIfPresent(coeffDict());
00195         C1_.readIfPresent(coeffDict());
00196         C2_.readIfPresent(coeffDict());
00197         sigmaEps_.readIfPresent(coeffDict());
00198 
00199         return true;
00200     }
00201     else
00202     {
00203         return false;
00204     }
00205 }
00206 
00207 
00208 void kEpsilon::correct()
00209 {
00210     RASModel::correct();
00211 
00212     if (!turbulence_)
00213     {
00214         return;
00215     }
00216 
00217     volScalarField G("RASModel::G", nut_*2*magSqr(symm(fvc::grad(U_))));
00218 
00219     // Update espsilon and G at the wall
00220     epsilon_.boundaryField().updateCoeffs();
00221 
00222     // Dissipation equation
00223     tmp<fvScalarMatrix> epsEqn
00224     (
00225         fvm::ddt(epsilon_)
00226       + fvm::div(phi_, epsilon_)
00227       - fvm::Sp(fvc::div(phi_), epsilon_)
00228       - fvm::laplacian(DepsilonEff(), epsilon_)
00229      ==
00230         C1_*G*epsilon_/k_
00231       - fvm::Sp(C2_*epsilon_/k_, epsilon_)
00232     );
00233 
00234     epsEqn().relax();
00235 
00236     epsEqn().boundaryManipulate(epsilon_.boundaryField());
00237 
00238     solve(epsEqn);
00239     bound(epsilon_, epsilon0_);
00240 
00241 
00242     // Turbulent kinetic energy equation
00243     tmp<fvScalarMatrix> kEqn
00244     (
00245         fvm::ddt(k_)
00246       + fvm::div(phi_, k_)
00247       - fvm::Sp(fvc::div(phi_), k_)
00248       - fvm::laplacian(DkEff(), k_)
00249      ==
00250         G
00251       - fvm::Sp(epsilon_/k_, k_)
00252     );
00253 
00254     kEqn().relax();
00255     solve(kEqn);
00256     bound(k_, k0_);
00257 
00258 
00259     // Re-calculate viscosity
00260     nut_ = Cmu_*sqr(k_)/epsilon_;
00261     nut_.correctBoundaryConditions();
00262 }
00263 
00264 
00265 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00266 
00267 } // End namespace RASModels
00268 } // End namespace incompressible
00269 } // End namespace Foam
00270 
00271 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines