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

qZeta.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 "qZeta.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(qZeta, 0);
00043 addToRunTimeSelectionTable(RASModel, qZeta, dictionary);
00044 
00045 // * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //
00046 
00047 tmp<volScalarField> qZeta::fMu() const
00048 {
00049     volScalarField Rt = q_*k_/(2.0*nu()*zeta_);
00050 
00051     if (anisotropic_)
00052     {
00053         return exp((-scalar(2.5) + Rt/20.0)/pow(scalar(1) + Rt/130.0, 3.0));
00054     }
00055     else
00056     {
00057         return
00058             exp(-6.0/sqr(scalar(1) + Rt/50.0))
00059            *(scalar(1) + 3.0*exp(-Rt/10.0));
00060     }
00061 }
00062 
00063 
00064 tmp<volScalarField> qZeta::f2() const
00065 {
00066     volScalarField Rt = q_*k_/(2.0*nu()*zeta_);
00067     return scalar(1) - 0.3*exp(-sqr(Rt));
00068 }
00069 
00070 
00071 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00072 
00073 qZeta::qZeta
00074 (
00075     const volVectorField& U,
00076     const surfaceScalarField& phi,
00077     transportModel& lamTransportModel
00078 )
00079 :
00080     RASModel(typeName, U, phi, lamTransportModel),
00081 
00082     Cmu_
00083     (
00084         dimensioned<scalar>::lookupOrAddToDict
00085         (
00086             "Cmu",
00087             coeffDict_,
00088             0.09
00089         )
00090     ),
00091     C1_
00092     (
00093         dimensioned<scalar>::lookupOrAddToDict
00094         (
00095             "C1",
00096             coeffDict_,
00097             1.44
00098         )
00099     ),
00100     C2_
00101     (
00102         dimensioned<scalar>::lookupOrAddToDict
00103         (
00104             "C2",
00105             coeffDict_,
00106             1.92
00107         )
00108     ),
00109     sigmaZeta_
00110     (
00111         dimensioned<scalar>::lookupOrAddToDict
00112         (
00113             "sigmaZeta",
00114             coeffDict_,
00115             1.3
00116         )
00117     ),
00118     anisotropic_
00119     (
00120         Switch::lookupOrAddToDict
00121         (
00122             "anisotropic",
00123             coeffDict_,
00124             false
00125         )
00126     ),
00127 
00128     k_
00129     (
00130         IOobject
00131         (
00132             "k",
00133             runTime_.timeName(),
00134             mesh_,
00135             IOobject::MUST_READ,
00136             IOobject::AUTO_WRITE
00137         ),
00138         mesh_
00139     ),
00140 
00141     epsilon_
00142     (
00143         IOobject
00144         (
00145             "epsilon",
00146             runTime_.timeName(),
00147             mesh_,
00148             IOobject::MUST_READ,
00149             IOobject::AUTO_WRITE
00150         ),
00151         mesh_
00152     ),
00153 
00154     q_
00155     (
00156         IOobject
00157         (
00158             "q",
00159             runTime_.timeName(),
00160             mesh_,
00161             IOobject::NO_READ,
00162             IOobject::AUTO_WRITE
00163         ),
00164         sqrt(k_),
00165         k_.boundaryField().types()
00166     ),
00167 
00168     zeta_
00169     (
00170         IOobject
00171         (
00172             "zeta",
00173             runTime_.timeName(),
00174             mesh_,
00175             IOobject::NO_READ,
00176             IOobject::AUTO_WRITE
00177         ),
00178         epsilon_/(2.0*q_),
00179         epsilon_.boundaryField().types()
00180     ),
00181 
00182     nut_
00183     (
00184         IOobject
00185         (
00186             "nut",
00187             runTime_.timeName(),
00188             mesh_,
00189             IOobject::NO_READ,
00190             IOobject::AUTO_WRITE
00191         ),
00192         autoCreateNut("nut", mesh_)
00193     )
00194 {
00195     nut_ = Cmu_*fMu()*sqr(k_)/(epsilon_ + epsilonSmall_);
00196     nut_.correctBoundaryConditions();
00197 
00198     printCoeffs();
00199 }
00200 
00201 
00202 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00203 
00204 tmp<volSymmTensorField> qZeta::R() const
00205 {
00206     return tmp<volSymmTensorField>
00207     (
00208         new volSymmTensorField
00209         (
00210             IOobject
00211             (
00212                 "R",
00213                 runTime_.timeName(),
00214                 mesh_,
00215                 IOobject::NO_READ,
00216                 IOobject::NO_WRITE
00217             ),
00218             ((2.0/3.0)*I)*k_ - nut_*twoSymm(fvc::grad(U_)),
00219             k_.boundaryField().types()
00220         )
00221     );
00222 }
00223 
00224 
00225 tmp<volSymmTensorField> qZeta::devReff() const
00226 {
00227     return tmp<volSymmTensorField>
00228     (
00229         new volSymmTensorField
00230         (
00231             IOobject
00232             (
00233                 "devRhoReff",
00234                 runTime_.timeName(),
00235                 mesh_,
00236                 IOobject::NO_READ,
00237                 IOobject::NO_WRITE
00238             ),
00239            -nuEff()*dev(twoSymm(fvc::grad(U_)))
00240         )
00241     );
00242 }
00243 
00244 
00245 tmp<fvVectorMatrix> qZeta::divDevReff(volVectorField& U) const
00246 {
00247     return
00248     (
00249       - fvm::laplacian(nuEff(), U)
00250       - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00251     );
00252 }
00253 
00254 
00255 bool qZeta::read()
00256 {
00257     if (RASModel::read())
00258     {
00259         Cmu_.readIfPresent(coeffDict());
00260         C1_.readIfPresent(coeffDict());
00261         C2_.readIfPresent(coeffDict());
00262         sigmaZeta_.readIfPresent(coeffDict());
00263         anisotropic_.readIfPresent("anisotropic", coeffDict());
00264 
00265         return true;
00266     }
00267     else
00268     {
00269         return false;
00270     }
00271 }
00272 
00273 
00274 void qZeta::correct()
00275 {
00276     RASModel::correct();
00277 
00278     if (!turbulence_)
00279     {
00280         return;
00281     }
00282 
00283     volScalarField S2 = 2*magSqr(symm(fvc::grad(U_)));
00284 
00285     volScalarField G("RASModel::G", nut_/(2.0*q_)*S2);
00286     volScalarField E = nu()*nut_/q_*fvc::magSqrGradGrad(U_);
00287 
00288 
00289     // Zeta equation
00290 
00291     tmp<fvScalarMatrix> zetaEqn
00292     (
00293         fvm::ddt(zeta_)
00294       + fvm::div(phi_, zeta_)
00295       - fvm::laplacian(DzetaEff(), zeta_)
00296      ==
00297         (2.0*C1_ - 1)*G*zeta_/q_
00298       - fvm::Sp((2.0*C2_ - dimensionedScalar(1.0))*f2()*zeta_/q_, zeta_)
00299       + E
00300     );
00301 
00302     zetaEqn().relax();
00303     solve(zetaEqn);
00304     bound(zeta_, epsilon0_/(2*sqrt(k0_)));
00305 
00306 
00307     // q equation
00308 
00309     tmp<fvScalarMatrix> qEqn
00310     (
00311         fvm::ddt(q_)
00312       + fvm::div(phi_, q_)
00313       - fvm::laplacian(DqEff(), q_)
00314      ==
00315         G - fvm::Sp(zeta_/q_, q_)
00316     );
00317 
00318     qEqn().relax();
00319     solve(qEqn);
00320     bound(q_, sqrt(k0_));
00321 
00322 
00323     // Re-calculate k and epsilon
00324     k_ = sqr(q_);
00325     k_.correctBoundaryConditions();
00326 
00327     epsilon_ = 2*q_*zeta_;
00328     epsilon_.correctBoundaryConditions();
00329 
00330 
00331     // Re-calculate viscosity
00332     nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
00333     nut_.correctBoundaryConditions();
00334 }
00335 
00336 
00337 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00338 
00339 } // End namespace RASModels
00340 } // End namespace incompressible
00341 } // End namespace Foam
00342 
00343 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines