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

LamBremhorstKE.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 "LamBremhorstKE.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(LamBremhorstKE, 0);
00043 addToRunTimeSelectionTable(RASModel, LamBremhorstKE, dictionary);
00044 
00045 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00046 
00047 LamBremhorstKE::LamBremhorstKE
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             "alphaEps",
00088             coeffDict_,
00089             1.3
00090         )
00091     ),
00092 
00093     k_
00094     (
00095         IOobject
00096         (
00097             "k",
00098             runTime_.timeName(),
00099             mesh_,
00100             IOobject::MUST_READ,
00101             IOobject::AUTO_WRITE
00102         ),
00103         mesh_
00104     ),
00105 
00106     epsilon_
00107     (
00108         IOobject
00109         (
00110             "epsilon",
00111             runTime_.timeName(),
00112             mesh_,
00113             IOobject::MUST_READ,
00114             IOobject::AUTO_WRITE
00115         ),
00116         mesh_
00117     ),
00118 
00119     y_(mesh_),
00120 
00121     Rt_(sqr(k_)/(nu()*epsilon_)),
00122 
00123     fMu_
00124     (
00125         sqr(scalar(1) - exp(-0.0165*(sqrt(k_)*y_/nu())))
00126        *(scalar(1) + 20.5/(Rt_ + SMALL))
00127     ),
00128 
00129     nut_
00130     (
00131         IOobject
00132         (
00133             "nut",
00134             runTime_.timeName(),
00135             mesh_,
00136             IOobject::NO_READ,
00137             IOobject::AUTO_WRITE
00138         ),
00139         autoCreateLowReNut("nut", mesh_)
00140     )
00141 {
00142     nut_ = Cmu_*fMu_*sqr(k_)/(epsilon_ + epsilonSmall_);
00143     nut_.correctBoundaryConditions();
00144 
00145     printCoeffs();
00146 }
00147 
00148 
00149 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00150 
00151 tmp<volSymmTensorField> LamBremhorstKE::R() const
00152 {
00153     return tmp<volSymmTensorField>
00154     (
00155         new volSymmTensorField
00156         (
00157             IOobject
00158             (
00159                 "R",
00160                 runTime_.timeName(),
00161                 mesh_,
00162                 IOobject::NO_READ,
00163                 IOobject::NO_WRITE
00164             ),
00165             ((2.0/3.0)*I)*k_ - nut_*twoSymm(fvc::grad(U_)),
00166             k_.boundaryField().types()
00167         )
00168     );
00169 }
00170 
00171 
00172 tmp<volSymmTensorField> LamBremhorstKE::devReff() const
00173 {
00174     return tmp<volSymmTensorField>
00175     (
00176         new volSymmTensorField
00177         (
00178             IOobject
00179             (
00180                 "devRhoReff",
00181                 runTime_.timeName(),
00182                 mesh_,
00183                 IOobject::NO_READ,
00184                 IOobject::NO_WRITE
00185             ),
00186            -nuEff()*dev(twoSymm(fvc::grad(U_)))
00187         )
00188     );
00189 }
00190 
00191 
00192 tmp<fvVectorMatrix> LamBremhorstKE::divDevReff(volVectorField& U) const
00193 {
00194     return
00195     (
00196       - fvm::laplacian(nuEff(), U)
00197       - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00198     );
00199 }
00200 
00201 
00202 bool LamBremhorstKE::read()
00203 {
00204     if (RASModel::read())
00205     {
00206         Cmu_.readIfPresent(coeffDict());
00207         C1_.readIfPresent(coeffDict());
00208         C2_.readIfPresent(coeffDict());
00209         sigmaEps_.readIfPresent(coeffDict());
00210 
00211         return true;
00212     }
00213     else
00214     {
00215         return false;
00216     }
00217 }
00218 
00219 
00220 void LamBremhorstKE::correct()
00221 {
00222     RASModel::correct();
00223 
00224     if (!turbulence_)
00225     {
00226         return;
00227     }
00228 
00229     if (mesh_.changing())
00230     {
00231         y_.correct();
00232     }
00233 
00234     volScalarField G("RASModel::G", nut_*2*magSqr(symm(fvc::grad(U_))));
00235 
00236 
00237     // Calculate parameters and coefficients for low-Reynolds number model
00238 
00239     Rt_ = sqr(k_)/(nu()*epsilon_);
00240     volScalarField Ry = sqrt(k_)*y_/nu();
00241 
00242     fMu_ = sqr(scalar(1) - exp(-0.0165*Ry))*(scalar(1) + 20.5/(Rt_ + SMALL));
00243 
00244     volScalarField f1 = scalar(1) + pow(0.05/(fMu_ + SMALL), 3);
00245     volScalarField f2 = scalar(1) - exp(-sqr(Rt_));
00246 
00247 
00248     // Dissipation equation
00249 
00250     tmp<fvScalarMatrix> epsEqn
00251     (
00252         fvm::ddt(epsilon_)
00253       + fvm::div(phi_, epsilon_)
00254       - fvm::laplacian(DepsilonEff(), epsilon_)
00255      ==
00256         C1_*f1*G*epsilon_/k_
00257       - fvm::Sp(C2_*f2*epsilon_/k_, epsilon_)
00258     );
00259 
00260     epsEqn().relax();
00261     solve(epsEqn);
00262     bound(epsilon_, epsilon0_);
00263 
00264 
00265     // Turbulent kinetic energy equation
00266 
00267     tmp<fvScalarMatrix> kEqn
00268     (
00269         fvm::ddt(k_)
00270       + fvm::div(phi_, k_)
00271       - fvm::laplacian(DkEff(), k_)
00272      ==
00273         G - fvm::Sp(epsilon_/k_, k_)
00274     );
00275 
00276     kEqn().relax();
00277     solve(kEqn);
00278     bound(k_, k0_);
00279 
00280 
00281     // Re-calculate viscosity
00282     nut_ == Cmu_*fMu_*sqr(k_)/epsilon_;
00283 }
00284 
00285 
00286 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00287 
00288 } // End namespace RASModels
00289 } // End namespace incompressible
00290 } // End namespace Foam
00291 
00292 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines