initial commit for version 1.5.x patch release
[OpenFOAM-1.5.x.git] / src / turbulenceModels / LES / compressible / lowReOneEqEddy / lowReOneEqEddy.C
bloba626068335caa63e9ca5df35ae5b70c24933ad91
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2008 OpenCFD Ltd.
6      \\/     M anipulation  |
7 -------------------------------------------------------------------------------
8 License
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 the
13     Free Software Foundation; either version 2 of the License, or (at your
14     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
19     for more details.
21     You should have received a copy of the GNU General Public License
22     along with OpenFOAM; if not, write to the Free Software Foundation,
23     Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301 USA
25 \*---------------------------------------------------------------------------*/
27 #include "lowReOneEqEddy.H"
28 #include "addToRunTimeSelectionTable.H"
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 namespace Foam
34 namespace compressible
36 namespace LESModels
39 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
41 defineTypeNameAndDebug(lowReOneEqEddy, 0);
42 addToRunTimeSelectionTable(LESModel, lowReOneEqEddy, dictionary);
44 // * * * * * * * * * * * * * * * * Constructors  * * * * * * * * * * * * * * //
46 // from components
47 lowReOneEqEddy::lowReOneEqEddy
49     const volScalarField& rho,
50     const volVectorField& U,
51     const surfaceScalarField& phi,
52     const basicThermo& thermoPhysicalModel
55     LESModel(typeName, rho, U, phi, thermoPhysicalModel),
56     GenEddyVisc(rho, U, phi, thermoPhysicalModel),
58     ck_
59     (
60         dimensioned<scalar>::lookupOrAddToDict
61         (
62             "ck",
63             coeffDict(),
64             0.07
65         )
66     ),
67     beta_
68     (
69         dimensioned<scalar>::lookupOrAddToDict
70         (
71             "beta",
72             coeffDict(),
73             0.01
74         )
75     )
77     printCoeffs();
81 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
83 void lowReOneEqEddy::correct(const tmp<volTensorField>& tgradU)
85     const volTensorField& gradU = tgradU();
87     GenEddyVisc::correct(gradU);
89     volScalarField divU = fvc::div(phi()/fvc::interpolate(rho()));
90     volScalarField G = 2*muSgs_*(gradU && dev(symm(gradU)));
92     solve
93     (
94         fvm::ddt(rho(), k_)
95       + fvm::div(phi(), k_)
96       - fvm::laplacian(DkEff(), k_)
97      ==
98         G
99       - fvm::SuSp(2.0/3.0*rho()*divU, k_)
100       - fvm::Sp(ce_*rho()*sqrt(k_)/delta(), k_)
101     );
103     bound(k_, k0());
105     // High Re eddy viscosity
106     muSgs_ = ck_*rho()*sqrt(k_)*delta();
108     // low Re no corrected eddy viscosity
109     muSgs_ -= (mu()/beta_)*(scalar(1) - exp(-beta_*muSgs_/mu()));
111     muSgs_.correctBoundaryConditions();
115 bool lowReOneEqEddy::read()
117     if (GenEddyVisc::read())
118     {
119         ck_.readIfPresent(coeffDict());
120         beta_.readIfPresent(coeffDict());
122         return true;
123     }
124     else
125     {
126         return false;
127     }
131 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
133 } // End namespace LESModels
134 } // End namespace compressible
135 } // End namespace Foam
137 // ************************************************************************* //