ENH: patchCloud: return pTraits<Type>::max for unfound points
[OpenFOAM-1.7.x.git] / src / turbulenceModels / compressible / LES / lowReOneEqEddy / lowReOneEqEddy.C
blob88e74e1a396165fc237d9ddbc1d352f60c758c7e
1 /*---------------------------------------------------------------------------*\
2   =========                 |
3   \\      /  F ield         | OpenFOAM: The Open Source CFD Toolbox
4    \\    /   O peration     |
5     \\  /    A nd           | Copyright (C) 1991-2010 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
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
19     for more details.
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 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 namespace Foam
33 namespace compressible
35 namespace LESModels
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
69     LESModel(typeName, rho, U, phi, thermoPhysicalModel),
70     GenEddyVisc(rho, U, phi, thermoPhysicalModel),
72     ck_
73     (
74         dimensioned<scalar>::lookupOrAddToDict
75         (
76             "ck",
77             coeffDict_,
78             0.07
79         )
80     ),
81     beta_
82     (
83         dimensioned<scalar>::lookupOrAddToDict
84         (
85             "beta",
86             coeffDict_,
87             0.01
88         )
89     )
91     updateSubGridScaleFields();
93     printCoeffs();
97 // * * * * * * * * * * * * * * * Member Functions  * * * * * * * * * * * * * //
99 void lowReOneEqEddy::correct(const tmp<volTensorField>& tgradU)
101     const volTensorField& gradU = tgradU();
103     GenEddyVisc::correct(gradU);
105     volScalarField divU = fvc::div(phi()/fvc::interpolate(rho()));
106     volScalarField G = 2*muSgs_*(gradU && dev(symm(gradU)));
108     solve
109     (
110         fvm::ddt(rho(), k_)
111       + fvm::div(phi(), k_)
112       - fvm::laplacian(DkEff(), k_)
113      ==
114         G
115       - fvm::SuSp(2.0/3.0*rho()*divU, k_)
116       - fvm::Sp(ce_*rho()*sqrt(k_)/delta(), k_)
117     );
119     bound(k_, k0());
121     updateSubGridScaleFields();
125 bool lowReOneEqEddy::read()
127     if (GenEddyVisc::read())
128     {
129         ck_.readIfPresent(coeffDict());
130         beta_.readIfPresent(coeffDict());
132         return true;
133     }
134     else
135     {
136         return false;
137     }
141 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
143 } // End namespace LESModels
144 } // End namespace compressible
145 } // End namespace Foam
147 // ************************************************************************* //