LienCubicKE.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 "LienCubicKE.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(LienCubicKE, 0);
43 addToRunTimeSelectionTable(RASModel, LienCubicKE, 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) + (2*kappa_/(pow(Cmu_, 0.75))/(yStar + small)));
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) + (2*kappa_/(pow(Cmu_, 0.75))/(yStar + small)))
71  );
72 
73  return
74  (Ceps2_*pow(Cmu_, 0.75))
75  *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
76 }
77 
78 
80 {
82 }
83 
84 
86 {
87  volSymmTensorField S(symm(gradU));
88  volTensorField W(skew(gradU));
89 
90  volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
91  volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
92 
93  volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
94  volScalarField fMu(this->fMu());
95 
96  nut_ = Cmu*fMu*sqr(k_)/epsilon_;
98 
100  fMu*k_
101  *(
102  // Quadratic terms
103  sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
104  *(
105  Cbeta1_*dev(innerSqr(S))
106  + Cbeta2_*twoSymm(S&W)
107  + Cbeta3_*dev(symm(W&W))
108  )
109 
110  // Cubic terms
111  - pow3(Cmu*k_/epsilon_)
112  *(
113  (Cgamma1_*magSqr(S) - Cgamma2_*magSqr(W))*S
114  + Cgamma4_*twoSymm((innerSqr(S)&W))
115  )
116  );
117 }
118 
119 
120 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
121 
123 (
124  const geometricOneField& alpha,
125  const geometricOneField& rho,
126  const volVectorField& U,
127  const surfaceScalarField& alphaRhoPhi,
128  const surfaceScalarField& phi,
129  const transportModel& transport,
130  const word& type
131 )
132 :
134  (
135  type,
136  alpha,
137  rho,
138  U,
139  alphaRhoPhi,
140  phi,
141  transport
142  ),
143 
144  Ceps1_
145  (
147  (
148  "Ceps1",
149  coeffDict_,
150  1.44
151  )
152  ),
153  Ceps2_
154  (
156  (
157  "Ceps2",
158  coeffDict_,
159  1.92
160  )
161  ),
162  sigmak_
163  (
165  (
166  "sigmak",
167  coeffDict_,
168  1.0
169  )
170  ),
171  sigmaEps_
172  (
174  (
175  "sigmaEps",
176  coeffDict_,
177  1.3
178  )
179  ),
180  Cmu1_
181  (
183  (
184  "Cmu1",
185  coeffDict_,
186  1.25
187  )
188  ),
189  Cmu2_
190  (
192  (
193  "Cmu2",
194  coeffDict_,
195  0.9
196  )
197  ),
198  Cbeta_
199  (
201  (
202  "Cbeta",
203  coeffDict_,
204  1000.0
205  )
206  ),
207  Cbeta1_
208  (
210  (
211  "Cbeta1",
212  coeffDict_,
213  3.0
214  )
215  ),
216  Cbeta2_
217  (
219  (
220  "Cbeta2",
221  coeffDict_,
222  15.0
223  )
224  ),
225  Cbeta3_
226  (
228  (
229  "Cbeta3",
230  coeffDict_,
231  -19.0
232  )
233  ),
234  Cgamma1_
235  (
237  (
238  "Cgamma1",
239  coeffDict_,
240  16.0
241  )
242  ),
243  Cgamma2_
244  (
246  (
247  "Cgamma2",
248  coeffDict_,
249  16.0
250  )
251  ),
252  Cgamma4_
253  (
255  (
256  "Cgamma4",
257  coeffDict_,
258  -80.0
259  )
260  ),
261  Cmu_
262  (
264  (
265  "Cmu",
266  coeffDict_,
267  0.09
268  )
269  ),
270  kappa_
271  (
273  (
274  "kappa",
275  coeffDict_,
276  0.41
277  )
278  ),
279  Anu_
280  (
282  (
283  "Anu",
284  coeffDict_,
285  0.0198
286  )
287  ),
288  AE_
289  (
291  (
292  "AE",
293  coeffDict_,
294  0.00375
295  )
296  ),
297 
298  k_
299  (
300  IOobject
301  (
302  IOobject::groupName("k", alphaRhoPhi.group()),
303  runTime_.timeName(),
304  mesh_,
307  ),
308  mesh_
309  ),
310 
311  epsilon_
312  (
313  IOobject
314  (
315  IOobject::groupName("epsilon", alphaRhoPhi.group()),
316  runTime_.timeName(),
317  mesh_,
320  ),
321  mesh_
322  ),
323 
324  y_(wallDist::New(mesh_).y())
325 {
326  bound(k_, kMin_);
327  bound(epsilon_, epsilonMin_);
328 
329  if (type == typeName)
330  {
331  printCoeffs(type);
332  }
333 }
334 
335 
336 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
337 
339 {
341  {
342  Ceps1_.readIfPresent(coeffDict());
343  Ceps2_.readIfPresent(coeffDict());
344  sigmak_.readIfPresent(coeffDict());
345  sigmaEps_.readIfPresent(coeffDict());
346  Cmu1_.readIfPresent(coeffDict());
347  Cmu2_.readIfPresent(coeffDict());
348  Cbeta_.readIfPresent(coeffDict());
349  Cbeta1_.readIfPresent(coeffDict());
350  Cbeta2_.readIfPresent(coeffDict());
351  Cbeta3_.readIfPresent(coeffDict());
352  Cgamma1_.readIfPresent(coeffDict());
353  Cgamma2_.readIfPresent(coeffDict());
354  Cgamma4_.readIfPresent(coeffDict());
355  Cmu_.readIfPresent(coeffDict());
356  kappa_.readIfPresent(coeffDict());
357  Anu_.readIfPresent(coeffDict());
358  AE_.readIfPresent(coeffDict());
359 
360  return true;
361  }
362  else
363  {
364  return false;
365  }
366 }
367 
368 
370 {
371  if (!turbulence_)
372  {
373  return;
374  }
375 
377 
378  tmp<volTensorField> tgradU = fvc::grad(U_);
379  const volTensorField& gradU = tgradU();
380 
382  (
383  GName(),
384  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
385  );
386 
387 
388  // Update epsilon and G at the wall
389  epsilon_.boundaryFieldRef().updateCoeffs();
390 
391  const volScalarField f2(this->f2());
392 
393  // Dissipation equation
394  tmp<fvScalarMatrix> epsEqn
395  (
397  + fvm::div(phi_, epsilon_)
399  ==
400  Ceps1_*G*epsilon_/k_
402  + E(f2)
403  );
404 
405  epsEqn.ref().relax();
406  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
407  solve(epsEqn);
408  bound(epsilon_, epsilonMin_);
409 
410 
411  // Turbulent kinetic energy equation
413  (
414  fvm::ddt(k_)
415  + fvm::div(phi_, k_)
416  - fvm::laplacian(DkEff(), k_)
417  ==
418  G
419  - fvm::Sp(epsilon_/k_, k_)
420  );
421 
422  kEqn.ref().relax();
423  solve(kEqn);
424  bound(k_, kMin_);
425 
426 
427  // Re-calculate viscosity and non-linear stress
428  correctNonlinearStress(gradU);
429 }
430 
431 
432 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
433 
434 } // End namespace RASModels
435 } // End namespace incompressible
436 } // End namespace Foam
437 
438 // ************************************************************************* //
const volScalarField & y_
Wall distance.
Definition: LienCubicKE.H:114
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensionedTensor skew(const dimensionedTensor &dt)
tmp< volScalarField > fMu() const
Definition: LienCubicKE.C:47
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...
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
virtual void correctNonlinearStress(const volTensorField &gradU)
Definition: LienCubicKE.C:85
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const GeometricField< Type, fvPatchField, volMesh > &)
static word groupName(Name name, const word &group)
LienCubicKE(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.
Definition: LienCubicKE.C:123
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual bool read()
Read RASProperties dictionary.
Definition: LienCubicKE.C:338
phi
Definition: correctPhi.H:3
RASModel< momentumTransportModel > RASModel
incompressible::RASModel ::transportModel transportModel
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
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.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LienCubicKE.H:168
U
Definition: pEqn.H:72
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LienCubicKE.C:369
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.
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
addToRunTimeSelectionTable(RASModel, kkLOmega, dictionary)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
const scalarList W(::W(thermo))
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
tmp< volScalarField > f2() const
Definition: LienCubicKE.C:57
void correctBoundaryConditions()
Correct boundary field.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: LienCubicKE.H:158
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
tmp< volScalarField > E(const volScalarField &f2) const
Definition: LienCubicKE.C:65
Namespace for OpenFOAM.