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-2023 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 
32 (
33  geometricOneField,
34  geometricOneField,
35  incompressibleMomentumTransportModel
36 )
37 
38 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
39 
40 namespace Foam
41 {
42 namespace incompressible
43 {
44 namespace RASModels
45 {
46 
47 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
48 
49 defineTypeNameAndDebug(LienCubicKE, 0);
51 (
52  RASincompressibleMomentumTransportModel,
53  LienCubicKE,
54  dictionary
55 );
56 
57 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
58 
59 tmp<volScalarField> LienCubicKE::boundEpsilon()
60 {
61  tmp<volScalarField> tCmuk2(Cmu_*sqr(k_));
62  epsilon_ = max(epsilon_, tCmuk2()/(nutMaxCoeff_*nu()));
63  return tCmuk2;
64 }
65 
66 
67 tmp<volScalarField> LienCubicKE::fMu() const
68 {
69  const volScalarField yStar(sqrt(k_)*y()/nu());
70 
71  return
72  (scalar(1) - exp(-Anu_*yStar))
73  *(scalar(1) + (2*kappa_/(pow(Cmu_, 0.75))/(yStar + small)));
74 }
75 
76 
77 tmp<volScalarField> LienCubicKE::f2() const
78 {
79  tmp<volScalarField> Rt = sqr(k_)/(nu()*epsilon_);
80 
81  return scalar(1) - 0.3*exp(-sqr(Rt));
82 }
83 
84 
85 tmp<volScalarField> LienCubicKE::E(const volScalarField& f2) const
86 {
87  const volScalarField yStar(sqrt(k_)*y()/nu());
88  const volScalarField le
89  (
90  kappa_*y()/(scalar(1) + (2*kappa_/(pow(Cmu_, 0.75))/(yStar + small)))
91  );
92 
93  return
94  (Ceps2_*pow(Cmu_, 0.75))
95  *(f2*sqrt(k_)*epsilon_/le)*exp(-AE_*sqr(yStar));
96 }
97 
98 
100 {
102 }
103 
104 
106 {
107  volSymmTensorField S(symm(gradU));
108  volTensorField W(skew(gradU));
109 
110  volScalarField sBar((k_/epsilon_)*sqrt(2.0)*mag(S));
111  volScalarField wBar((k_/epsilon_)*sqrt(2.0)*mag(W));
112 
113  volScalarField Cmu((2.0/3.0)/(Cmu1_ + sBar + Cmu2_*wBar));
114  volScalarField fMu(this->fMu());
115 
116  nut_ = fMu*boundEpsilon()/epsilon_;
118 
120  fMu*k_
121  *(
122  // Quadratic terms
123  sqr(k_/epsilon_)/(Cbeta_ + pow3(sBar))
124  *(
126  + Cbeta2_*twoSymm(S&W)
127  + Cbeta3_*dev(symm(W&W))
128  )
129 
130  // Cubic terms
131  - pow3(Cmu*k_/epsilon_)
132  *(
134  + Cgamma4_*twoSymm((innerSqr(S)&W))
135  )
136  );
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
141 
143 (
144  const geometricOneField& alpha,
145  const geometricOneField& rho,
146  const volVectorField& U,
147  const surfaceScalarField& alphaRhoPhi,
148  const surfaceScalarField& phi,
149  const viscosity& viscosity,
150  const word& type
151 )
152 :
153  nonlinearEddyViscosity<incompressible::RASModel>
154  (
155  type,
156  alpha,
157  rho,
158  U,
159  alphaRhoPhi,
160  phi,
161  viscosity
162  ),
163 
164  Ceps1_
165  (
167  (
168  "Ceps1",
169  coeffDict_,
170  1.44
171  )
172  ),
173  Ceps2_
174  (
176  (
177  "Ceps2",
178  coeffDict_,
179  1.92
180  )
181  ),
182  sigmak_
183  (
185  (
186  "sigmak",
187  coeffDict_,
188  1.0
189  )
190  ),
191  sigmaEps_
192  (
194  (
195  "sigmaEps",
196  coeffDict_,
197  1.3
198  )
199  ),
200  Cmu1_
201  (
203  (
204  "Cmu1",
205  coeffDict_,
206  1.25
207  )
208  ),
209  Cmu2_
210  (
212  (
213  "Cmu2",
214  coeffDict_,
215  0.9
216  )
217  ),
218  Cbeta_
219  (
221  (
222  "Cbeta",
223  coeffDict_,
224  1000.0
225  )
226  ),
227  Cbeta1_
228  (
230  (
231  "Cbeta1",
232  coeffDict_,
233  3.0
234  )
235  ),
236  Cbeta2_
237  (
239  (
240  "Cbeta2",
241  coeffDict_,
242  15.0
243  )
244  ),
245  Cbeta3_
246  (
248  (
249  "Cbeta3",
250  coeffDict_,
251  -19.0
252  )
253  ),
254  Cgamma1_
255  (
257  (
258  "Cgamma1",
259  coeffDict_,
260  16.0
261  )
262  ),
263  Cgamma2_
264  (
266  (
267  "Cgamma2",
268  coeffDict_,
269  16.0
270  )
271  ),
272  Cgamma4_
273  (
275  (
276  "Cgamma4",
277  coeffDict_,
278  -80.0
279  )
280  ),
281  Cmu_
282  (
284  (
285  "Cmu",
286  coeffDict_,
287  0.09
288  )
289  ),
290  kappa_
291  (
293  (
294  "kappa",
295  coeffDict_,
296  0.41
297  )
298  ),
299  Anu_
300  (
302  (
303  "Anu",
304  coeffDict_,
305  0.0198
306  )
307  ),
308  AE_
309  (
311  (
312  "AE",
313  coeffDict_,
314  0.00375
315  )
316  ),
317 
318  k_
319  (
320  IOobject
321  (
322  this->groupName("k"),
323  runTime_.name(),
324  mesh_,
327  ),
328  mesh_
329  ),
330 
331  epsilon_
332  (
333  IOobject
334  (
335  this->groupName("epsilon"),
336  runTime_.name(),
337  mesh_,
340  ),
341  mesh_
342  )
343 {
344  bound(k_, kMin_);
345  boundEpsilon();
346 
347  if (type == typeName)
348  {
349  printCoeffs(type);
350  }
351 }
352 
353 
354 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
355 
356 bool LienCubicKE::read()
357 {
359  {
360  Ceps1_.readIfPresent(coeffDict());
361  Ceps2_.readIfPresent(coeffDict());
362  sigmak_.readIfPresent(coeffDict());
363  sigmaEps_.readIfPresent(coeffDict());
364  Cmu1_.readIfPresent(coeffDict());
365  Cmu2_.readIfPresent(coeffDict());
366  Cbeta_.readIfPresent(coeffDict());
367  Cbeta1_.readIfPresent(coeffDict());
368  Cbeta2_.readIfPresent(coeffDict());
369  Cbeta3_.readIfPresent(coeffDict());
370  Cgamma1_.readIfPresent(coeffDict());
371  Cgamma2_.readIfPresent(coeffDict());
372  Cgamma4_.readIfPresent(coeffDict());
373  Cmu_.readIfPresent(coeffDict());
374  kappa_.readIfPresent(coeffDict());
375  Anu_.readIfPresent(coeffDict());
376  AE_.readIfPresent(coeffDict());
377 
378  return true;
379  }
380  else
381  {
382  return false;
383  }
384 }
385 
386 
388 {
389  if (!turbulence_)
390  {
391  return;
392  }
393 
395 
396  tmp<volTensorField> tgradU = fvc::grad(U_);
397  const volTensorField& gradU = tgradU();
398 
400  (
401  GName(),
402  (nut_*twoSymm(gradU) - nonlinearStress_) && gradU
403  );
404 
405 
406  // Update epsilon and G at the wall
408 
409  const volScalarField f2(this->f2());
410 
411  // Dissipation equation
412  tmp<fvScalarMatrix> epsEqn
413  (
415  + fvm::div(phi_, epsilon_)
417  ==
420  + E(f2)
421  );
422 
423  epsEqn.ref().relax();
424  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
425  solve(epsEqn);
426  boundEpsilon();
427 
428 
429  // Turbulent kinetic energy equation
430  tmp<fvScalarMatrix> kEqn
431  (
432  fvm::ddt(k_)
433  + fvm::div(phi_, k_)
434  - fvm::laplacian(DkEff(), k_)
435  ==
436  G
437  - fvm::Sp(epsilon_/k_, k_)
438  );
439 
440  kEqn.ref().relax();
441  solve(kEqn);
442  bound(k_, kMin_);
443 
444 
445  // Re-calculate viscosity and non-linear stress
446  correctNonlinearStress(gradU);
447 }
448 
449 
450 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
451 
452 } // End namespace RASModels
453 } // End namespace incompressible
454 } // End namespace Foam
455 
456 // ************************************************************************* //
scalar y
makeMomentumTransportModelTypes(geometricOneField, geometricOneField, incompressibleMomentumTransportModel)
Definition: LienCubicKE.C:32
Bound the given scalar field where it is below the specified minimum.
void updateCoeffs()
Update the boundary condition coefficients.
void correctBoundaryConditions()
Correct boundary field.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
static dimensioned< scalar > lookupOrAddToDict(const word &, dictionary &, const dimensionSet &dims=dimless, const scalar &defaultValue=pTraits< scalar >::zero)
Construct from dictionary, with default value.
bool readIfPresent(const dictionary &)
Update the value of dimensioned<Type> if found in the dictionary.
virtual bool read()=0
Re-read model coefficients if they have changed.
Definition: eddyViscosity.C:74
virtual void correct()=0
Solve the turbulence equations and correct the turbulence viscosity.
LienCubicKE(const geometricOneField &alpha, const geometricOneField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual void correctNut()
Correct the eddy-viscosity nut.
tmp< volScalarField > DkEff() const
Return the effective diffusivity for k.
Definition: LienCubicKE.H:158
tmp< volScalarField > f2() const
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
tmp< volScalarField > E(const volScalarField &f2) const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
tmp< volScalarField > DepsilonEff() const
Return the effective diffusivity for epsilon.
Definition: LienCubicKE.H:168
virtual void correctNonlinearStress(const volTensorField &gradU)
virtual bool read()
Read RASProperties dictionary.
tmp< volScalarField > fMu() const
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< fvMatrix< Type > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:65
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
SurfaceField< scalar > surfaceScalarField
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:68
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
scalarList W(const fluidMulticomponentThermo &thermo)
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
dimensionedSymmTensor innerSqr(const dimensionedSymmTensor &dt)
VolField< symmTensor > volSymmTensorField
Definition: volFieldsFwd.H:67
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
Info<< "Reading field U\n"<< endl;volVectorField U(IOobject("U", runTime.name(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Creating face flux\n"<< endl;surfaceScalarField phi(IOobject("phi", runTime.name(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), mesh, dimensionedScalar(mesh.Sf().dimensions() *U.dimensions(), 0));autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))