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

LienCubicKE.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 "LienCubicKE.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(LienCubicKE, 0);
00043 addToRunTimeSelectionTable(RASModel, LienCubicKE, dictionary);
00044 
00045 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00046 
00047 LienCubicKE::LienCubicKE
00048 (
00049     const volVectorField& U,
00050     const surfaceScalarField& phi,
00051     transportModel& lamTransportModel
00052 )
00053 :
00054     RASModel(typeName, U, phi, lamTransportModel),
00055 
00056     C1_
00057     (
00058         dimensioned<scalar>::lookupOrAddToDict
00059         (
00060             "C1",
00061             coeffDict_,
00062             1.44
00063         )
00064     ),
00065     C2_
00066     (
00067         dimensioned<scalar>::lookupOrAddToDict
00068         (
00069             "C2",
00070             coeffDict_,
00071             1.92
00072         )
00073     ),
00074     sigmak_
00075     (
00076         dimensioned<scalar>::lookupOrAddToDict
00077         (
00078             "sigmak",
00079             coeffDict_,
00080             1.0
00081         )
00082     ),
00083     sigmaEps_
00084     (
00085         dimensioned<scalar>::lookupOrAddToDict
00086         (
00087             "sigmaEps",
00088             coeffDict_,
00089             1.3
00090         )
00091     ),
00092     A1_
00093     (
00094         dimensioned<scalar>::lookupOrAddToDict
00095         (
00096             "A1",
00097             coeffDict_,
00098             1.25
00099         )
00100     ),
00101     A2_
00102     (
00103         dimensioned<scalar>::lookupOrAddToDict
00104         (
00105             "A2",
00106             coeffDict_,
00107             1000.0
00108         )
00109     ),
00110     Ctau1_
00111     (
00112         dimensioned<scalar>::lookupOrAddToDict
00113         (
00114             "Ctau1",
00115             coeffDict_,
00116             -4.0
00117         )
00118     ),
00119     Ctau2_
00120     (
00121         dimensioned<scalar>::lookupOrAddToDict
00122         (
00123             "Ctau2",
00124             coeffDict_,
00125             13.0
00126         )
00127     ),
00128     Ctau3_
00129     (
00130         dimensioned<scalar>::lookupOrAddToDict
00131         (
00132             "Ctau3",
00133             coeffDict_,
00134             -2.0
00135         )
00136     ),
00137     alphaKsi_
00138     (
00139         dimensioned<scalar>::lookupOrAddToDict
00140         (
00141             "alphaKsi",
00142             coeffDict_,
00143             0.9
00144         )
00145     ),
00146 
00147     k_
00148     (
00149         IOobject
00150         (
00151             "k",
00152             runTime_.timeName(),
00153             mesh_,
00154             IOobject::NO_READ,
00155             IOobject::AUTO_WRITE
00156         ),
00157         autoCreateK("k", mesh_)
00158     ),
00159     epsilon_
00160     (
00161         IOobject
00162         (
00163             "epsilon",
00164             runTime_.timeName(),
00165             mesh_,
00166             IOobject::NO_READ,
00167             IOobject::AUTO_WRITE
00168         ),
00169         autoCreateEpsilon("epsilon", mesh_)
00170     ),
00171 
00172     gradU_(fvc::grad(U)),
00173     eta_(k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ + gradU_.T())))),
00174     ksi_(k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ - gradU_.T())))),
00175     Cmu_(2.0/(3.0*(A1_ + eta_ + alphaKsi_*ksi_))),
00176     fEta_(A2_ + pow(eta_, 3.0)),
00177 
00178     C5viscosity_
00179     (
00180       - 2.0*pow3(Cmu_)*pow4(k_)/pow3(epsilon_)
00181        *(
00182             magSqr(gradU_ + gradU_.T())
00183           - magSqr(gradU_ - gradU_.T())
00184         )
00185     ),
00186 
00187     nut_
00188     (
00189         IOobject
00190         (
00191             "nut",
00192             runTime_.timeName(),
00193             mesh_,
00194             IOobject::NO_READ,
00195             IOobject::AUTO_WRITE
00196         ),
00197         autoCreateNut("nut", mesh_)
00198     ),
00199 
00200     nonlinearStress_
00201     (
00202         "nonlinearStress",
00203         // quadratic terms
00204         symm
00205         (
00206             pow(k_, 3.0)/sqr(epsilon_)
00207            *(
00208                 Ctau1_/fEta_
00209                *(
00210                     (gradU_ & gradU_)
00211                   + (gradU_ & gradU_)().T()
00212                 )
00213               + Ctau2_/fEta_*(gradU_ & gradU_.T())
00214               + Ctau3_/fEta_*(gradU_.T() & gradU_)
00215             )
00216             // cubic term C4
00217           - 20.0*pow(k_, 4.0)/pow(epsilon_, 3.0)
00218            *pow(Cmu_, 3.0)
00219            *(
00220                 ((gradU_ & gradU_) & gradU_.T())
00221               + ((gradU_ & gradU_.T()) & gradU_.T())
00222               - ((gradU_.T() & gradU_) & gradU_)
00223               - ((gradU_.T() & gradU_.T()) & gradU_)
00224             )
00225         )
00226     )
00227 {
00228     nut_ = Cmu_*sqr(k_)/(epsilon_ + epsilonSmall_) + C5viscosity_;
00229     nut_.correctBoundaryConditions();
00230 
00231     printCoeffs();
00232 }
00233 
00234 
00235 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00236 
00237 tmp<volSymmTensorField> LienCubicKE::R() const
00238 {
00239     return tmp<volSymmTensorField>
00240     (
00241         new volSymmTensorField
00242         (
00243             IOobject
00244             (
00245                 "R",
00246                 runTime_.timeName(),
00247                 mesh_,
00248                 IOobject::NO_READ,
00249                 IOobject::NO_WRITE
00250             ),
00251             ((2.0/3.0)*I)*k_ - nut_*twoSymm(gradU_) + nonlinearStress_,
00252             k_.boundaryField().types()
00253         )
00254     );
00255 }
00256 
00257 
00258 tmp<volSymmTensorField> LienCubicKE::devReff() const
00259 {
00260     return tmp<volSymmTensorField>
00261     (
00262         new volSymmTensorField
00263         (
00264             IOobject
00265             (
00266                 "devRhoReff",
00267                 runTime_.timeName(),
00268                 mesh_,
00269                 IOobject::NO_READ,
00270                 IOobject::NO_WRITE
00271             ),
00272            -nuEff()*dev(twoSymm(fvc::grad(U_))) + nonlinearStress_
00273         )
00274     );
00275 }
00276 
00277 
00278 tmp<fvVectorMatrix> LienCubicKE::divDevReff(volVectorField& U) const
00279 {
00280     return
00281     (
00282         fvc::div(nonlinearStress_)
00283       - fvm::laplacian(nuEff(), U)
00284       - fvc::div(nuEff()*dev(fvc::grad(U)().T()))
00285     );
00286 }
00287 
00288 
00289 bool LienCubicKE::read()
00290 {
00291     if (RASModel::read())
00292     {
00293         C1_.readIfPresent(coeffDict());
00294         C2_.readIfPresent(coeffDict());
00295         sigmak_.readIfPresent(coeffDict());
00296         sigmaEps_.readIfPresent(coeffDict());
00297         A1_.readIfPresent(coeffDict());
00298         A2_.readIfPresent(coeffDict());
00299         Ctau1_.readIfPresent(coeffDict());
00300         Ctau2_.readIfPresent(coeffDict());
00301         Ctau3_.readIfPresent(coeffDict());
00302         alphaKsi_.readIfPresent(coeffDict());
00303 
00304         return true;
00305     }
00306     else
00307     {
00308         return false;
00309     }
00310 }
00311 
00312 
00313 void LienCubicKE::correct()
00314 {
00315     RASModel::correct();
00316 
00317     if (!turbulence_)
00318     {
00319         return;
00320     }
00321 
00322     gradU_ = fvc::grad(U_);
00323 
00324     // generation term
00325     volScalarField S2 = symm(gradU_) && gradU_;
00326 
00327     volScalarField G
00328     (
00329         "RASModel::G",
00330         Cmu_*sqr(k_)/epsilon_*S2 - (nonlinearStress_ && gradU_)
00331     );
00332 
00333     // Update espsilon and G at the wall
00334     epsilon_.boundaryField().updateCoeffs();
00335 
00336     // Dissipation equation
00337     tmp<fvScalarMatrix> epsEqn
00338     (
00339         fvm::ddt(epsilon_)
00340       + fvm::div(phi_, epsilon_)
00341       - fvm::laplacian(DepsilonEff(), epsilon_)
00342       ==
00343         C1_*G*epsilon_/k_
00344       - fvm::Sp(C2_*epsilon_/k_, epsilon_)
00345     );
00346 
00347     epsEqn().relax();
00348 
00349     epsEqn().boundaryManipulate(epsilon_.boundaryField());
00350 
00351     solve(epsEqn);
00352     bound(epsilon_, epsilon0_);
00353 
00354 
00355     // Turbulent kinetic energy equation
00356 
00357     tmp<fvScalarMatrix> kEqn
00358     (
00359         fvm::ddt(k_)
00360       + fvm::div(phi_, k_)
00361       - fvm::laplacian(DkEff(), k_)
00362       ==
00363         G
00364       - fvm::Sp(epsilon_/k_, k_)
00365     );
00366 
00367     kEqn().relax();
00368     solve(kEqn);
00369     bound(k_, k0_);
00370 
00371 
00372     // Re-calculate viscosity
00373 
00374     eta_ = k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ + gradU_.T())));
00375     ksi_ = k_/epsilon_*sqrt(2.0*magSqr(0.5*(gradU_ - gradU_.T())));
00376     Cmu_ = 2.0/(3.0*(A1_ + eta_ + alphaKsi_*ksi_));
00377     fEta_ = A2_ + pow(eta_, 3.0);
00378 
00379     C5viscosity_ =
00380         - 2.0*pow(Cmu_, 3.0)*pow(k_, 4.0)/pow(epsilon_, 3.0)
00381        *(magSqr(gradU_ + gradU_.T()) - magSqr(gradU_ - gradU_.T()));
00382 
00383     nut_ = Cmu_*sqr(k_)/epsilon_ + C5viscosity_;
00384     nut_.correctBoundaryConditions();
00385 
00386     nonlinearStress_ = symm
00387     (
00388         // quadratic terms
00389         pow(k_, 3.0)/sqr(epsilon_)*
00390         (
00391             Ctau1_/fEta_*
00392             (
00393                 (gradU_ & gradU_)
00394               + (gradU_ & gradU_)().T()
00395             )
00396           + Ctau2_/fEta_*(gradU_ & gradU_.T())
00397           + Ctau3_/fEta_*(gradU_.T() & gradU_)
00398         )
00399         // cubic term C4
00400       - 20.0*pow(k_, 4.0)/pow(epsilon_, 3.0)
00401        *pow(Cmu_, 3.0)
00402        *(
00403             ((gradU_ & gradU_) & gradU_.T())
00404           + ((gradU_ & gradU_.T()) & gradU_.T())
00405           - ((gradU_.T() & gradU_) & gradU_)
00406           - ((gradU_.T() & gradU_.T()) & gradU_)
00407         )
00408     );
00409 }
00410 
00411 
00412 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00413 
00414 } // End namespace RASModels
00415 } // End namespace incompressible
00416 } // End namespace Foam
00417 
00418 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines