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

kOmegaSST.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 "kOmegaSST.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(kOmegaSST, 0);
00043 addToRunTimeSelectionTable(RASModel, kOmegaSST, dictionary);
00044 
00045 // * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * * //
00046 
00047 tmp<volScalarField> kOmegaSST::F1(const volScalarField& CDkOmega) const
00048 {
00049     volScalarField CDkOmegaPlus = max
00050     (
00051         CDkOmega,
00052         dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
00053     );
00054 
00055     volScalarField arg1 = min
00056     (
00057         min
00058         (
00059             max
00060             (
00061                 (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
00062                 scalar(500)*(mu()/rho_)/(sqr(y_)*omega_)
00063             ),
00064             (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
00065         ),
00066         scalar(10)
00067     );
00068 
00069     return tanh(pow4(arg1));
00070 }
00071 
00072 tmp<volScalarField> kOmegaSST::F2() const
00073 {
00074     volScalarField arg2 = min
00075     (
00076         max
00077         (
00078             (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
00079             scalar(500)*(mu()/rho_)/(sqr(y_)*omega_)
00080         ),
00081         scalar(100)
00082     );
00083 
00084     return tanh(sqr(arg2));
00085 }
00086 
00087 
00088 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00089 
00090 kOmegaSST::kOmegaSST
00091 (
00092     const volScalarField& rho,
00093     const volVectorField& U,
00094     const surfaceScalarField& phi,
00095     const basicThermo& thermophysicalModel
00096 )
00097 :
00098     RASModel(typeName, rho, U, phi, thermophysicalModel),
00099 
00100     alphaK1_
00101     (
00102         dimensioned<scalar>::lookupOrAddToDict
00103         (
00104             "alphaK1",
00105             coeffDict_,
00106             0.85034
00107         )
00108     ),
00109     alphaK2_
00110     (
00111         dimensioned<scalar>::lookupOrAddToDict
00112         (
00113             "alphaK2",
00114             coeffDict_,
00115             1.0
00116         )
00117     ),
00118     alphaOmega1_
00119     (
00120         dimensioned<scalar>::lookupOrAddToDict
00121         (
00122             "alphaOmega1",
00123             coeffDict_,
00124             0.5
00125         )
00126     ),
00127     alphaOmega2_
00128     (
00129         dimensioned<scalar>::lookupOrAddToDict
00130         (
00131             "alphaOmega2",
00132             coeffDict_,
00133             0.85616
00134         )
00135     ),
00136     Prt_
00137     (
00138         dimensioned<scalar>::lookupOrAddToDict
00139         (
00140             "Prt",
00141             coeffDict_,
00142             1.0
00143         )
00144     ),
00145     gamma1_
00146     (
00147         dimensioned<scalar>::lookupOrAddToDict
00148         (
00149             "gamma1",
00150             coeffDict_,
00151             0.5532
00152         )
00153     ),
00154     gamma2_
00155     (
00156         dimensioned<scalar>::lookupOrAddToDict
00157         (
00158             "gamma2",
00159             coeffDict_,
00160             0.4403
00161         )
00162     ),
00163     beta1_
00164     (
00165         dimensioned<scalar>::lookupOrAddToDict
00166         (
00167             "beta1",
00168             coeffDict_,
00169             0.075
00170         )
00171     ),
00172     beta2_
00173     (
00174         dimensioned<scalar>::lookupOrAddToDict
00175         (
00176             "beta2",
00177             coeffDict_,
00178             0.0828
00179         )
00180     ),
00181     betaStar_
00182     (
00183         dimensioned<scalar>::lookupOrAddToDict
00184         (
00185             "betaStar",
00186             coeffDict_,
00187             0.09
00188         )
00189     ),
00190     a1_
00191     (
00192         dimensioned<scalar>::lookupOrAddToDict
00193         (
00194             "a1",
00195             coeffDict_,
00196             0.31
00197         )
00198     ),
00199     c1_
00200     (
00201         dimensioned<scalar>::lookupOrAddToDict
00202         (
00203             "c1",
00204             coeffDict_,
00205             10.0
00206         )
00207     ),
00208 
00209     y_(mesh_),
00210 
00211     k_
00212     (
00213         IOobject
00214         (
00215             "k",
00216             runTime_.timeName(),
00217             mesh_,
00218             IOobject::NO_READ,
00219             IOobject::AUTO_WRITE
00220         ),
00221         autoCreateK("k", mesh_)
00222     ),
00223     omega_
00224     (
00225         IOobject
00226         (
00227             "omega",
00228             runTime_.timeName(),
00229             mesh_,
00230             IOobject::NO_READ,
00231             IOobject::AUTO_WRITE
00232         ),
00233         autoCreateOmega("omega", mesh_)
00234     ),
00235     mut_
00236     (
00237         IOobject
00238         (
00239             "mut",
00240             runTime_.timeName(),
00241             mesh_,
00242             IOobject::NO_READ,
00243             IOobject::AUTO_WRITE
00244         ),
00245         autoCreateMut("mut", mesh_)
00246     ),
00247     alphat_
00248     (
00249         IOobject
00250         (
00251             "alphat",
00252             runTime_.timeName(),
00253             mesh_,
00254             IOobject::NO_READ,
00255             IOobject::AUTO_WRITE
00256         ),
00257         autoCreateAlphat("alphat", mesh_)
00258     )
00259 {
00260     bound(omega_, omega0_);
00261 
00262     mut_ =
00263     (
00264         a1_*rho_*k_
00265       / max
00266         (
00267             a1_*omega_,
00268             F2()*sqrt(2.0*magSqr(symm(fvc::grad(U_))))
00269         )
00270     );
00271     mut_.correctBoundaryConditions();
00272 
00273     alphat_ = mut_/Prt_;
00274     alphat_.correctBoundaryConditions();
00275 
00276     printCoeffs();
00277 }
00278 
00279 
00280 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00281 
00282 tmp<volSymmTensorField> kOmegaSST::R() const
00283 {
00284     return tmp<volSymmTensorField>
00285     (
00286         new volSymmTensorField
00287         (
00288             IOobject
00289             (
00290                 "R",
00291                 runTime_.timeName(),
00292                 mesh_,
00293                 IOobject::NO_READ,
00294                 IOobject::NO_WRITE
00295             ),
00296             ((2.0/3.0)*I)*k_ - (mut_/rho_)*dev(twoSymm(fvc::grad(U_))),
00297             k_.boundaryField().types()
00298         )
00299     );
00300 }
00301 
00302 
00303 tmp<volSymmTensorField> kOmegaSST::devRhoReff() const
00304 {
00305     return tmp<volSymmTensorField>
00306     (
00307         new volSymmTensorField
00308         (
00309             IOobject
00310             (
00311                 "devRhoReff",
00312                 runTime_.timeName(),
00313                 mesh_,
00314                 IOobject::NO_READ,
00315                 IOobject::NO_WRITE
00316             ),
00317             -muEff()*dev(twoSymm(fvc::grad(U_)))
00318         )
00319     );
00320 }
00321 
00322 
00323 tmp<fvVectorMatrix> kOmegaSST::divDevRhoReff(volVectorField& U) const
00324 {
00325     return
00326     (
00327       - fvm::laplacian(muEff(), U) - fvc::div(muEff()*dev2(fvc::grad(U)().T()))
00328     );
00329 }
00330 
00331 
00332 bool kOmegaSST::read()
00333 {
00334     if (RASModel::read())
00335     {
00336         alphaK1_.readIfPresent(coeffDict());
00337         alphaK2_.readIfPresent(coeffDict());
00338         alphaOmega1_.readIfPresent(coeffDict());
00339         alphaOmega2_.readIfPresent(coeffDict());
00340         Prt_.readIfPresent(coeffDict());
00341         gamma1_.readIfPresent(coeffDict());
00342         gamma2_.readIfPresent(coeffDict());
00343         beta1_.readIfPresent(coeffDict());
00344         beta2_.readIfPresent(coeffDict());
00345         betaStar_.readIfPresent(coeffDict());
00346         a1_.readIfPresent(coeffDict());
00347         c1_.readIfPresent(coeffDict());
00348 
00349         return true;
00350     }
00351     else
00352     {
00353         return false;
00354     }
00355 }
00356 
00357 
00358 void kOmegaSST::correct()
00359 {
00360     if (!turbulence_)
00361     {
00362         // Re-calculate viscosity
00363         mut_ =
00364             a1_*rho_*k_
00365            /max(a1_*omega_, F2()*sqrt(2.0*magSqr(symm(fvc::grad(U_)))));
00366         mut_.correctBoundaryConditions();
00367 
00368         // Re-calculate thermal diffusivity
00369         alphat_ = mut_/Prt_;
00370         alphat_.correctBoundaryConditions();
00371 
00372         return;
00373     }
00374 
00375     RASModel::correct();
00376 
00377     volScalarField divU = fvc::div(phi_/fvc::interpolate(rho_));
00378 
00379     if (mesh_.changing())
00380     {
00381         y_.correct();
00382     }
00383 
00384     if (mesh_.moving())
00385     {
00386         divU += fvc::div(mesh_.phi());
00387     }
00388 
00389     tmp<volTensorField> tgradU = fvc::grad(U_);
00390     volScalarField S2 = magSqr(symm(tgradU()));
00391     volScalarField GbyMu = (tgradU() && dev(twoSymm(tgradU())));
00392     volScalarField G("RASModel::G", mut_*GbyMu);
00393     tgradU.clear();
00394 
00395     // Update omega and G at the wall
00396     omega_.boundaryField().updateCoeffs();
00397 
00398     volScalarField CDkOmega =
00399         (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_;
00400 
00401     volScalarField F1 = this->F1(CDkOmega);
00402     volScalarField rhoGammaF1 = rho_*gamma(F1);
00403 
00404     // Turbulent frequency equation
00405     tmp<fvScalarMatrix> omegaEqn
00406     (
00407         fvm::ddt(rho_, omega_)
00408       + fvm::div(phi_, omega_)
00409       - fvm::laplacian(DomegaEff(F1), omega_)
00410      ==
00411         rhoGammaF1*GbyMu
00412       - fvm::SuSp((2.0/3.0)*rhoGammaF1*divU, omega_)
00413       - fvm::Sp(rho_*beta(F1)*omega_, omega_)
00414       - fvm::SuSp
00415         (
00416             rho_*(F1 - scalar(1))*CDkOmega/omega_,
00417             omega_
00418         )
00419     );
00420 
00421     omegaEqn().relax();
00422 
00423     omegaEqn().boundaryManipulate(omega_.boundaryField());
00424 
00425     solve(omegaEqn);
00426     bound(omega_, omega0_);
00427 
00428     // Turbulent kinetic energy equation
00429     tmp<fvScalarMatrix> kEqn
00430     (
00431         fvm::ddt(rho_, k_)
00432       + fvm::div(phi_, k_)
00433       - fvm::laplacian(DkEff(F1), k_)
00434      ==
00435         min(G, (c1_*betaStar_)*rho_*k_*omega_)
00436       - fvm::SuSp(2.0/3.0*rho_*divU, k_)
00437       - fvm::Sp(rho_*betaStar_*omega_, k_)
00438     );
00439 
00440     kEqn().relax();
00441     solve(kEqn);
00442     bound(k_, k0_);
00443 
00444 
00445     // Re-calculate viscosity
00446     mut_ = a1_*rho_*k_/max(a1_*omega_, F2()*sqrt(2.0*S2));
00447     mut_.correctBoundaryConditions();
00448 
00449     // Re-calculate thermal diffusivity
00450     alphat_ = mut_/Prt_;
00451     alphat_.correctBoundaryConditions();
00452 }
00453 
00454 
00455 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00456 
00457 } // End namespace RASModels
00458 } // End namespace compressible
00459 } // End namespace Foam
00460 
00461 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines