LienLeschziner.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 "LienLeschziner.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(LienLeschziner, 0);
43 addToRunTimeSelectionTable(RASModel, LienLeschziner, dictionary);
44 
45 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
46 
48 {
49  const volScalarField yStar(sqrt(k_)*y_/nu());
50 
51  return
52  (scalar(1) - exp(-Anu_*yStar))
53  /((scalar(1) + small) - exp(-Aeps_*yStar));
54 }
55 
56 
58 {
59  tmp<volScalarField> Rt = sqr(k_)/(nu()*epsilon_);
60 
61  return scalar(1) - 0.3*exp(-sqr(Rt));
62 }
63 
64 
66 {
67  const volScalarField yStar(sqrt(k_)*y_/nu());
68  const volScalarField le
69  (
70  kappa_*y_*((scalar(1) + small) - exp(-Aeps_*yStar))
71  );
72 
73  return
74  (Ceps2_*pow(Cmu_, 0.75))
75  *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
76 }
77 
78 
80 {
81  nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
83 }
84 
85 
86 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
87 
89 (
90  const geometricOneField& alpha,
91  const geometricOneField& rho,
92  const volVectorField& U,
93  const surfaceScalarField& alphaRhoPhi,
94  const surfaceScalarField& phi,
95  const transportModel& transport,
96  const word& type
97 )
98 :
100  (
101  type,
102  alpha,
103  rho,
104  U,
105  alphaRhoPhi,
106  phi,
107  transport
108  ),
109 
110  Ceps1_
111  (
113  (
114  "Ceps1",
115  coeffDict_,
116  1.44
117  )
118  ),
119  Ceps2_
120  (
122  (
123  "Ceps2",
124  coeffDict_,
125  1.92
126  )
127  ),
128  sigmak_
129  (
131  (
132  "sigmak",
133  coeffDict_,
134  1.0
135  )
136  ),
137  sigmaEps_
138  (
140  (
141  "sigmaEps",
142  coeffDict_,
143  1.3
144  )
145  ),
146  Cmu_
147  (
149  (
150  "Cmu",
151  coeffDict_,
152  0.09
153  )
154  ),
155  kappa_
156  (
158  (
159  "kappa",
160  coeffDict_,
161  0.41
162  )
163  ),
164  Anu_
165  (
167  (
168  "Anu",
169  coeffDict_,
170  0.016
171  )
172  ),
173  Aeps_
174  (
176  (
177  "Aeps",
178  coeffDict_,
179  0.263
180  )
181  ),
182  AE_
183  (
185  (
186  "AE",
187  coeffDict_,
188  0.00222
189  )
190  ),
191 
192  k_
193  (
194  IOobject
195  (
196  IOobject::groupName("k", alphaRhoPhi.group()),
197  runTime_.timeName(),
198  mesh_,
201  ),
202  mesh_
203  ),
204 
205  epsilon_
206  (
207  IOobject
208  (
209  IOobject::groupName("epsilon", alphaRhoPhi.group()),
210  runTime_.timeName(),
211  mesh_,
214  ),
215  mesh_
216  ),
217 
218  y_(wallDist::New(mesh_).y())
219 {
220  bound(k_, kMin_);
221  bound(epsilon_, epsilonMin_);
222 
223  if (type == typeName)
224  {
225  printCoeffs(type);
226  }
227 }
228 
229 
230 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
231 
233 {
235  {
236  Ceps1_.readIfPresent(coeffDict());
237  Ceps2_.readIfPresent(coeffDict());
238  sigmak_.readIfPresent(coeffDict());
239  sigmaEps_.readIfPresent(coeffDict());
240  Cmu_.readIfPresent(coeffDict());
241  kappa_.readIfPresent(coeffDict());
242  Anu_.readIfPresent(coeffDict());
243  Aeps_.readIfPresent(coeffDict());
244  AE_.readIfPresent(coeffDict());
245 
246  return true;
247  }
248  else
249  {
250  return false;
251  }
252 }
253 
254 
256 {
257  if (!turbulence_)
258  {
259  return;
260  }
261 
263 
264  tmp<volTensorField> tgradU = fvc::grad(U_);
266  (
267  GName(),
268  nut_*(tgradU() && twoSymm(tgradU()))
269  );
270  tgradU.clear();
271 
272  // Update epsilon and G at the wall
273  epsilon_.boundaryFieldRef().updateCoeffs();
274 
275  const volScalarField f2(this->f2());
276 
277  // Dissipation equation
278  tmp<fvScalarMatrix> epsEqn
279  (
281  + fvm::div(phi_, epsilon_)
283  ==
284  Ceps1_*G*epsilon_/k_
286  + E(f2)
287  );
288 
289  epsEqn.ref().relax();
290  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
291  solve(epsEqn);
292  bound(epsilon_, epsilonMin_);
293 
294 
295  // Turbulent kinetic energy equation
297  (
298  fvm::ddt(k_)
299  + fvm::div(phi_, k_)
300  - fvm::laplacian(DkEff(), k_)
301  ==
302  G
303  - fvm::Sp(epsilon_/k_, k_)
304  );
305 
306  kEqn.ref().relax();
307  solve(kEqn);
308  bound(k_, kMin_);
309 
310  correctNut();
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 } // End namespace RASModels
317 } // End namespace incompressible
318 } // End namespace Foam
319 
320 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
tmp< volScalarField > E(const volScalarField &f2) const
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
Generic dimensioned Type class.
const dimensionedScalar G
Newtonian constant of gravitation.
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)
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)
incompressible::RASModel ::transportModel transportModel
Definition: eddyViscosity.H:72
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual bool read()
Read RASProperties dictionary.
phi
Definition: correctPhi.H:3
RASModel< momentumTransportModel > RASModel
const volScalarField & y_
Wall distance.
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.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
U
Definition: pEqn.H:72
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
void correctBoundaryConditions()
Correct boundary field.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
LienLeschziner(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Namespace for OpenFOAM.