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 | Copyright (C) 2011-2016 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
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.
15 
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.
20 
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/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "LamBremhorstKE.H"
27 #include "wallDist.H"
28 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace incompressible
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
41 
42 defineTypeNameAndDebug(LamBremhorstKE, 0);
43 addToRunTimeSelectionTable(RASModel, LamBremhorstKE, dictionary);
44 
45 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
46 
47 tmp<volScalarField> LamBremhorstKE::Rt() const
48 {
49  return sqr(k_)/(nu()*epsilon_);
50 }
51 
52 
53 tmp<volScalarField> LamBremhorstKE::fMu(const volScalarField& Rt) const
54 {
55  tmp<volScalarField> Ry(sqrt(k_)*y_/nu());
56  return sqr(scalar(1) - exp(-0.0165*Ry))*(scalar(1) + 20.5/(Rt + SMALL));
57 }
58 
59 
60 tmp<volScalarField> LamBremhorstKE::f1(const volScalarField& fMu) const
61 {
62  return scalar(1) + pow3(0.05/(fMu + SMALL));
63 }
64 
65 
66 tmp<volScalarField> LamBremhorstKE::f2(const volScalarField& Rt) const
67 {
68  return scalar(1) - exp(-sqr(Rt));
69 }
70 
71 
73 {
74  nut_ = Cmu_*fMu*sqr(k_)/epsilon_;
76 }
77 
78 
79 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
80 
82 {
83  correctNut(fMu(Rt()));
84 }
85 
86 
87 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
88 
89 LamBremhorstKE::LamBremhorstKE
90 (
91  const geometricOneField& alpha,
92  const geometricOneField& rho,
93  const volVectorField& U,
94  const surfaceScalarField& alphaRhoPhi,
95  const surfaceScalarField& phi,
96  const transportModel& transport,
97  const word& propertiesName,
98  const word& type
99 )
100 :
102  (
103  type,
104  alpha,
105  rho,
106  U,
107  alphaRhoPhi,
108  phi,
109  transport,
110  propertiesName
111  ),
112 
113  Cmu_
114  (
116  (
117  "Cmu",
118  coeffDict_,
119  0.09
120  )
121  ),
122  Ceps1_
123  (
125  (
126  "Ceps1",
127  coeffDict_,
128  1.44
129  )
130  ),
131  Ceps2_
132  (
134  (
135  "Ceps2",
136  coeffDict_,
137  1.92
138  )
139  ),
140  sigmaEps_
141  (
143  (
144  "alphaEps",
145  coeffDict_,
146  1.3
147  )
148  ),
149 
150  k_
151  (
152  IOobject
153  (
154  IOobject::groupName("k", U.group()),
155  runTime_.timeName(),
156  mesh_,
159  ),
160  mesh_
161  ),
162 
163  epsilon_
164  (
165  IOobject
166  (
167  IOobject::groupName("epsilon", U.group()),
168  runTime_.timeName(),
169  mesh_,
172  ),
173  mesh_
174  ),
175 
176  y_(wallDist::New(mesh_).y())
177 {
178  bound(k_, kMin_);
179  bound(epsilon_, epsilonMin_);
180 
181  if (type == typeName)
182  {
183  printCoeffs(type);
184  }
185 }
186 
187 
188 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
189 
191 {
193  {
194  Cmu_.readIfPresent(coeffDict());
195  Ceps1_.readIfPresent(coeffDict());
196  Ceps2_.readIfPresent(coeffDict());
197  sigmaEps_.readIfPresent(coeffDict());
198 
199  return true;
200  }
201  else
202  {
203  return false;
204  }
205 }
206 
207 
209 {
210  if (!turbulence_)
211  {
212  return;
213  }
214 
216 
217  tmp<volTensorField> tgradU = fvc::grad(U_);
218  volScalarField G(GName(), nut_*(twoSymm(tgradU()) && tgradU()));
219  tgradU.clear();
220 
221  // Update epsilon and G at the wall
222  epsilon_.boundaryFieldRef().updateCoeffs();
223 
224  const volScalarField Rt(this->Rt());
225  const volScalarField fMu(this->fMu(Rt));
226 
227  // Dissipation equation
228  tmp<fvScalarMatrix> epsEqn
229  (
231  + fvm::div(phi_, epsilon_)
233  ==
234  Ceps1_*f1(fMu)*G*epsilon_/k_
235  - fvm::Sp(Ceps2_*f2(Rt)*epsilon_/k_, epsilon_)
236  );
237 
238  epsEqn.ref().relax();
239  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
240  solve(epsEqn);
241  bound(epsilon_, epsilonMin_);
242 
243  // Turbulent kinetic energy equation
245  (
246  fvm::ddt(k_)
247  + fvm::div(phi_, k_)
248  - fvm::laplacian(DkEff(), k_)
249  ==
250  G - fvm::Sp(epsilon_/k_, k_)
251  );
252 
253  kEqn.ref().relax();
254  solve(kEqn);
255  bound(k_, kMin_);
256 
257  correctNut(fMu);
258 }
259 
260 
261 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
262 
263 } // End namespace RASModels
264 } // End namespace incompressible
265 } // End namespace Foam
266 
267 // ************************************************************************* //
surfaceScalarField & phi
virtual bool read()
Read RASProperties dictionary.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
U
Definition: pEqn.H:83
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
const dimensionedScalar G
Newtonian constant of gravitation.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
RASModel< turbulenceModel > RASModel
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Macros for easy insertion into run-time selection tables.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:75
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
scalar y
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Bound the given scalar field if it has gone unbounded.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
const volScalarField & y_
Wall distance.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
tensor Ry(const scalar &omega)
Rotational transformation tensor about the y-axis by omega radians.
Definition: transform.H:95
dimensionedScalar pow3(const dimensionedScalar &ds)
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
void correctBoundaryConditions()
Correct boundary field.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
A class for managing temporary objects.
Definition: PtrList.H:53
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
volScalarField & nu
Namespace for OpenFOAM.