1 /*---------------------------------------------------------------------------*\
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
5 \\ / A nd | Copyright (C) 2011 OpenFOAM Foundation
7 -------------------------------------------------------------------------------
9 This file is part of OpenFOAM.
11 OpenFOAM is free software: you can redistribute it and/or modify it
12 under the terms of the GNU General Public License as published by
13 the Free Software Foundation, either version 3 of the License, or
14 (at your option) any later version.
16 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 You should have received a copy of the GNU General Public License
22 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
24 \*---------------------------------------------------------------------------*/
26 #include "lowReOneEqEddy.H"
27 #include "addToRunTimeSelectionTable.H"
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 namespace compressible
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
40 defineTypeNameAndDebug(lowReOneEqEddy, 0);
41 addToRunTimeSelectionTable(LESModel, lowReOneEqEddy, dictionary);
43 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
45 void lowReOneEqEddy::updateSubGridScaleFields()
47 // High Re eddy viscosity
48 muSgs_ = ck_*rho()*sqrt(k_)*delta();
50 // low Re no corrected eddy viscosity
51 muSgs_ -= (mu()/beta_)*(scalar(1) - exp(-beta_*muSgs_/mu()));
52 muSgs_.correctBoundaryConditions();
54 alphaSgs_ = muSgs_/Prt_;
55 alphaSgs_.correctBoundaryConditions();
59 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
61 lowReOneEqEddy::lowReOneEqEddy
63 const volScalarField& rho,
64 const volVectorField& U,
65 const surfaceScalarField& phi,
66 const basicThermo& thermoPhysicalModel,
67 const word& turbulenceModelName,
71 LESModel(modelName, rho, U, phi, thermoPhysicalModel, turbulenceModelName),
72 GenEddyVisc(rho, U, phi, thermoPhysicalModel),
89 dimensioned<scalar>::lookupOrAddToDict
98 dimensioned<scalar>::lookupOrAddToDict
106 updateSubGridScaleFields();
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
114 void lowReOneEqEddy::correct(const tmp<volTensorField>& tgradU)
116 const volTensorField& gradU = tgradU();
118 GenEddyVisc::correct(gradU);
120 volScalarField divU(fvc::div(phi()/fvc::interpolate(rho())));
121 volScalarField G(2*muSgs_*(gradU && dev(symm(gradU))));
123 tmp<fvScalarMatrix> kEqn
126 + fvm::div(phi(), k_)
127 - fvm::laplacian(DkEff(), k_)
130 - fvm::SuSp(2.0/3.0*rho()*divU, k_)
131 - fvm::Sp(ce_*rho()*sqrt(k_)/delta(), k_)
139 updateSubGridScaleFields();
143 bool lowReOneEqEddy::read()
145 if (GenEddyVisc::read())
147 ck_.readIfPresent(coeffDict());
148 beta_.readIfPresent(coeffDict());
159 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
161 } // End namespace LESModels
162 } // End namespace compressible
163 } // End namespace Foam
165 // ************************************************************************* //