LienCubicKE.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 "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 {
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& propertiesName,
131  const word& type
132 )
133 :
135  (
136  type,
137  alpha,
138  rho,
139  U,
140  alphaRhoPhi,
141  phi,
142  transport,
143  propertiesName
144  ),
145 
146  Ceps1_
147  (
149  (
150  "Ceps1",
151  coeffDict_,
152  1.44
153  )
154  ),
155  Ceps2_
156  (
158  (
159  "Ceps2",
160  coeffDict_,
161  1.92
162  )
163  ),
164  sigmak_
165  (
167  (
168  "sigmak",
169  coeffDict_,
170  1.0
171  )
172  ),
173  sigmaEps_
174  (
176  (
177  "sigmaEps",
178  coeffDict_,
179  1.3
180  )
181  ),
182  Cmu1_
183  (
185  (
186  "Cmu1",
187  coeffDict_,
188  1.25
189  )
190  ),
191  Cmu2_
192  (
194  (
195  "Cmu2",
196  coeffDict_,
197  0.9
198  )
199  ),
200  Cbeta_
201  (
203  (
204  "Cbeta",
205  coeffDict_,
206  1000.0
207  )
208  ),
209  Cbeta1_
210  (
212  (
213  "Cbeta1",
214  coeffDict_,
215  3.0
216  )
217  ),
218  Cbeta2_
219  (
221  (
222  "Cbeta2",
223  coeffDict_,
224  15.0
225  )
226  ),
227  Cbeta3_
228  (
230  (
231  "Cbeta3",
232  coeffDict_,
233  -19.0
234  )
235  ),
236  Cgamma1_
237  (
239  (
240  "Cgamma1",
241  coeffDict_,
242  16.0
243  )
244  ),
245  Cgamma2_
246  (
248  (
249  "Cgamma2",
250  coeffDict_,
251  16.0
252  )
253  ),
254  Cgamma4_
255  (
257  (
258  "Cgamma4",
259  coeffDict_,
260  -80.0
261  )
262  ),
263  Cmu_
264  (
266  (
267  "Cmu",
268  coeffDict_,
269  0.09
270  )
271  ),
272  kappa_
273  (
275  (
276  "kappa",
277  coeffDict_,
278  0.41
279  )
280  ),
281  Anu_
282  (
284  (
285  "Anu",
286  coeffDict_,
287  0.0198
288  )
289  ),
290  AE_
291  (
293  (
294  "AE",
295  coeffDict_,
296  0.00375
297  )
298  ),
299 
300  k_
301  (
302  IOobject
303  (
304  IOobject::groupName("k", U.group()),
305  runTime_.timeName(),
306  mesh_,
309  ),
310  mesh_
311  ),
312 
313  epsilon_
314  (
315  IOobject
316  (
317  IOobject::groupName("epsilon", U.group()),
318  runTime_.timeName(),
319  mesh_,
322  ),
323  mesh_
324  ),
325 
326  y_(wallDist::New(mesh_).y())
327 {
328  bound(k_, kMin_);
329  bound(epsilon_, epsilonMin_);
330 
331  if (type == typeName)
332  {
333  printCoeffs(type);
334  }
335 }
336 
337 
338 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
339 
341 {
343  {
344  Ceps1_.readIfPresent(coeffDict());
345  Ceps2_.readIfPresent(coeffDict());
346  sigmak_.readIfPresent(coeffDict());
347  sigmaEps_.readIfPresent(coeffDict());
348  Cmu1_.readIfPresent(coeffDict());
349  Cmu2_.readIfPresent(coeffDict());
350  Cbeta_.readIfPresent(coeffDict());
351  Cbeta1_.readIfPresent(coeffDict());
352  Cbeta2_.readIfPresent(coeffDict());
353  Cbeta3_.readIfPresent(coeffDict());
354  Cgamma1_.readIfPresent(coeffDict());
355  Cgamma2_.readIfPresent(coeffDict());
356  Cgamma4_.readIfPresent(coeffDict());
357  Cmu_.readIfPresent(coeffDict());
358  kappa_.readIfPresent(coeffDict());
359  Anu_.readIfPresent(coeffDict());
360  AE_.readIfPresent(coeffDict());
361 
362  return true;
363  }
364  else
365  {
366  return false;
367  }
368 }
369 
370 
372 {
373  if (!turbulence_)
374  {
375  return;
376  }
377 
379 
380  tmp<volTensorField> tgradU = fvc::grad(U_);
381  const volTensorField& gradU = tgradU();
382 
384  (
385  GName(),
386  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
387  );
388 
389 
390  // Update epsilon and G at the wall
391  epsilon_.boundaryFieldRef().updateCoeffs();
392 
393  const volScalarField f2(this->f2());
394 
395  // Dissipation equation
396  tmp<fvScalarMatrix> epsEqn
397  (
399  + fvm::div(phi_, epsilon_)
401  ==
402  Ceps1_*G*epsilon_/k_
404  + E(f2)
405  );
406 
407  epsEqn.ref().relax();
408  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
409  solve(epsEqn);
410  bound(epsilon_, epsilonMin_);
411 
412 
413  // Turbulent kinetic energy equation
415  (
416  fvm::ddt(k_)
417  + fvm::div(phi_, k_)
418  - fvm::laplacian(DkEff(), k_)
419  ==
420  G
421  - fvm::Sp(epsilon_/k_, k_)
422  );
423 
424  kEqn.ref().relax();
425  solve(kEqn);
426  bound(k_, kMin_);
427 
428 
429  // Re-calculate viscosity and non-linear stress
430  correctNonlinearStress(gradU);
431 }
432 
433 
434 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
435 
436 } // End namespace RASModels
437 } // End namespace incompressible
438 } // End namespace Foam
439 
440 // ************************************************************************* //
const volScalarField & y_
Wall distance.
Definition: LienCubicKE.H:117
surfaceScalarField & phi
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
dimensionedTensor skew(const dimensionedTensor &dt)
tmp< volScalarField > fMu() const
Definition: LienCubicKE.C:47
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
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
LienCubicKE(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: LienCubicKE.C:123
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
virtual bool read()
Read RASProperties dictionary.
Definition: LienCubicKE.C:340
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.
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
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:171
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LienCubicKE.C:371
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
tmp< volScalarField > f2() const
Definition: LienCubicKE.C:57
void correctBoundaryConditions()
Correct boundary field.
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:162
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
tmp< volScalarField > E(const volScalarField &f2) const
Definition: LienCubicKE.C:65
Namespace for OpenFOAM.