LamBremhorstKE.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2011-2017 OpenFOAM Foundation
9 Copyright (C) 2019-2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
29#include "LamBremhorstKE.H"
30#include "wallDist.H"
31#include "bound.H"
33
34// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace incompressible
39{
40namespace RASModels
41{
42
43// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
44
47
48// * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
49
50tmp<volScalarField> LamBremhorstKE::Rt() const
51{
52 return sqr(k_)/(nu()*epsilon_);
53}
54
55
57{
59 return sqr(scalar(1) - exp(-0.0165*Ry))*(scalar(1) + 20.5/(Rt + SMALL));
60}
61
62
63tmp<volScalarField> LamBremhorstKE::f1(const volScalarField& fMu) const
64{
65 return scalar(1) + pow3(0.05/(fMu + SMALL));
66}
67
68
69tmp<volScalarField> LamBremhorstKE::f2(const volScalarField& Rt) const
70{
71 return scalar(1) - exp(-sqr(Rt));
72}
73
74
76{
77 nut_ = Cmu_*fMu*sqr(k_)/epsilon_;
79}
80
81
82// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
83
85{
86 correctNut(fMu(Rt()));
87}
88
89
90// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
91
93(
96 const volVectorField& U,
97 const surfaceScalarField& alphaRhoPhi,
99 const transportModel& transport,
100 const word& propertiesName,
101 const word& type
102)
103:
104 eddyViscosity<incompressible::RASModel>
105 (
106 type,
107 alpha,
108 rho,
109 U,
110 alphaRhoPhi,
111 phi,
112 transport,
113 propertiesName
114 ),
115
116 Cmu_
117 (
118 dimensioned<scalar>::getOrAddToDict
119 (
120 "Cmu",
121 coeffDict_,
122 0.09
123 )
124 ),
125 Ceps1_
126 (
127 dimensioned<scalar>::getOrAddToDict
128 (
129 "Ceps1",
130 coeffDict_,
131 1.44
132 )
133 ),
134 Ceps2_
135 (
136 dimensioned<scalar>::getOrAddToDict
137 (
138 "Ceps2",
139 coeffDict_,
140 1.92
141 )
142 ),
143 sigmaEps_
144 (
145 dimensioned<scalar>::getOrAddToDict
146 (
147 "alphaEps",
148 coeffDict_,
149 1.3
150 )
151 ),
152
153 k_
154 (
156 (
157 IOobject::groupName("k", alphaRhoPhi.group()),
158 runTime_.timeName(),
159 mesh_,
160 IOobject::MUST_READ,
161 IOobject::AUTO_WRITE
162 ),
163 mesh_
164 ),
165
166 epsilon_
167 (
169 (
170 IOobject::groupName("epsilon", alphaRhoPhi.group()),
171 runTime_.timeName(),
172 mesh_,
173 IOobject::MUST_READ,
174 IOobject::AUTO_WRITE
175 ),
176 mesh_
177 ),
178
179 y_(wallDist::New(mesh_).y())
180{
181 bound(k_, kMin_);
182 bound(epsilon_, epsilonMin_);
183
184 if (type == typeName)
185 {
186 printCoeffs(type);
187 }
188}
189
190
191// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192
194{
196 {
197 Cmu_.readIfPresent(coeffDict());
198 Ceps1_.readIfPresent(coeffDict());
199 Ceps2_.readIfPresent(coeffDict());
200 sigmaEps_.readIfPresent(coeffDict());
201
202 return true;
203 }
204
205 return false;
206}
207
208
210{
211 if (!turbulence_)
212 {
213 return;
214 }
215
217
218 tmp<volTensorField> tgradU = fvc::grad(U_);
219 volScalarField G(GName(), nut_*(twoSymm(tgradU()) && tgradU()));
220 tgradU.clear();
221
222 // Update epsilon and G at the wall
224
225 const volScalarField Rt(this->Rt());
226 const volScalarField fMu(this->fMu(Rt));
227
228 // Dissipation equation
230 (
232 + fvm::div(phi_, epsilon_)
234 ==
235 Ceps1_*f1(fMu)*G*epsilon_/k_
237 );
238
239 epsEqn.ref().relax();
240 epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
241 solve(epsEqn);
242 bound(epsilon_, epsilonMin_);
243
244 // Turbulent kinetic energy equation
246 (
247 fvm::ddt(k_)
248 + fvm::div(phi_, k_)
250 ==
251 G - fvm::Sp(epsilon_/k_, k_)
252 );
253
254 kEqn.ref().relax();
255 solve(kEqn);
256 bound(k_, kMin_);
257
258 // Update nut with latest available k,epsilon
259 correctNut();
260}
261
262
263// * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
264
265} // End namespace RASModels
266} // End namespace incompressible
267} // End namespace Foam
268
269// ************************************************************************* //
scalar y
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Bound the given scalar field if it has gone unbounded.
surfaceScalarField & phi
void updateCoeffs()
Update the boundary condition coefficients.
Boundary & boundaryFieldRef(const bool updateAccessTime=true)
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Generic dimensioned Type class.
bool readIfPresent(const dictionary &dict)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:58
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
Lam and Bremhorst low-Reynolds number k-epsilon turbulence model for incompressible flows.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
const volScalarField & y_
Wall distance.
virtual bool read()
Re-read model coefficients if they have changed.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
tmp< volScalarField > f2() const
Definition: LienCubicKE.C:60
tmp< volScalarField > fMu() const
Definition: LienCubicKE.C:50
BasicTurbulenceModel::transportModel transportModel
A class for managing temporary objects.
Definition: tmp.H:65
void clear() const noexcept
Definition: tmpI.H:287
T & ref() const
Definition: tmpI.H:227
Interface to run-time selectable methods to calculate the distance-to-wall and normal-to-wall fields.
Definition: wallDist.H:78
A class for handling words, derived from Foam::string.
Definition: word.H:68
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
U
Definition: pEqn.H:72
word timeName
Definition: getTimeIndex.H:3
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
autoPtr< BasicCompressibleTurbulenceModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const typename BasicCompressibleTurbulenceModel::transportModel &transport, const word &propertiesName)
RASModel< turbulenceModel > RASModel
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:82
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
tensor Ry(const scalar omega)
Rotational transformation tensor about the y-axis by omega radians.
Definition: transform.H:99
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:598
dimensionedScalar sqrt(const dimensionedScalar &ds)
volScalarField & nu
volScalarField & alpha
CEqn solve()