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

lowReOneEqEddy.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 "lowReOneEqEddy.H"
00027 #include <OpenFOAM/addToRunTimeSelectionTable.H>
00028 
00029 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00030 
00031 namespace Foam
00032 {
00033 namespace compressible
00034 {
00035 namespace LESModels
00036 {
00037 
00038 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00039 
00040 defineTypeNameAndDebug(lowReOneEqEddy, 0);
00041 addToRunTimeSelectionTable(LESModel, lowReOneEqEddy, dictionary);
00042 
00043 // * * * * * * * * * * * * * Private Member Functions  * * * * * * * * * * * //
00044 
00045 void lowReOneEqEddy::updateSubGridScaleFields()
00046 {
00047     // High Re eddy viscosity
00048     muSgs_ = ck_*rho()*sqrt(k_)*delta();
00049 
00050     // low Re no corrected eddy viscosity
00051     muSgs_ -= (mu()/beta_)*(scalar(1) - exp(-beta_*muSgs_/mu()));
00052     muSgs_.correctBoundaryConditions();
00053 
00054     alphaSgs_ = muSgs_/Prt_;
00055     alphaSgs_.correctBoundaryConditions();
00056 }
00057 
00058 
00059 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
00060 
00061 lowReOneEqEddy::lowReOneEqEddy
00062 (
00063     const volScalarField& rho,
00064     const volVectorField& U,
00065     const surfaceScalarField& phi,
00066     const basicThermo& thermoPhysicalModel
00067 )
00068 :
00069     LESModel(typeName, rho, U, phi, thermoPhysicalModel),
00070     GenEddyVisc(rho, U, phi, thermoPhysicalModel),
00071 
00072     ck_
00073     (
00074         dimensioned<scalar>::lookupOrAddToDict
00075         (
00076             "ck",
00077             coeffDict_,
00078             0.07
00079         )
00080     ),
00081     beta_
00082     (
00083         dimensioned<scalar>::lookupOrAddToDict
00084         (
00085             "beta",
00086             coeffDict_,
00087             0.01
00088         )
00089     )
00090 {
00091     updateSubGridScaleFields();
00092 
00093     printCoeffs();
00094 }
00095 
00096 
00097 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
00098 
00099 void lowReOneEqEddy::correct(const tmp<volTensorField>& tgradU)
00100 {
00101     const volTensorField& gradU = tgradU();
00102 
00103     GenEddyVisc::correct(gradU);
00104 
00105     volScalarField divU = fvc::div(phi()/fvc::interpolate(rho()));
00106     volScalarField G = 2*muSgs_*(gradU && dev(symm(gradU)));
00107 
00108     solve
00109     (
00110         fvm::ddt(rho(), k_)
00111       + fvm::div(phi(), k_)
00112       - fvm::laplacian(DkEff(), k_)
00113      ==
00114         G
00115       - fvm::SuSp(2.0/3.0*rho()*divU, k_)
00116       - fvm::Sp(ce_*rho()*sqrt(k_)/delta(), k_)
00117     );
00118 
00119     bound(k_, k0());
00120 
00121     updateSubGridScaleFields();
00122 }
00123 
00124 
00125 bool lowReOneEqEddy::read()
00126 {
00127     if (GenEddyVisc::read())
00128     {
00129         ck_.readIfPresent(coeffDict());
00130         beta_.readIfPresent(coeffDict());
00131 
00132         return true;
00133     }
00134     else
00135     {
00136         return false;
00137     }
00138 }
00139 
00140 
00141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
00142 
00143 } // End namespace LESModels
00144 } // End namespace compressible
00145 } // End namespace Foam
00146 
00147 // ************************ vim: set sw=4 sts=4 et: ************************ //
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines