kkLOmega.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-2024 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 "kkLOmega.H"
27 #include "bound.H"
28 #include "wallDist.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(kkLOmega, 0);
51 (
52  RASincompressibleMomentumTransportModel,
53  kkLOmega,
54  dictionary
55 );
56 
57 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
58 
59 tmp<volScalarField> kkLOmega::fv(const volScalarField& Ret) const
60 {
61  return(1 - exp(-sqrt(Ret)/Av_));
62 }
63 
64 
65 tmp<volScalarField> kkLOmega::fINT() const
66 {
67  return
68  (
69  min
70  (
71  kt_/(Cint_*(kl_ + kt_ + kMin_)),
73  )
74  );
75 }
76 
77 
78 tmp<volScalarField> kkLOmega::fSS(const volScalarField& Omega) const
79 {
80  return(exp(-sqr(Css_*nu()*Omega/(kt_ + kMin_))));
81 }
82 
83 
84 tmp<volScalarField> kkLOmega::Cmu(const volScalarField& S) const
85 {
86  return(1/(A0_ + As_*(S/(omega_ + omegaMin_))));
87 }
88 
89 
90 tmp<volScalarField> kkLOmega::BetaTS(const volScalarField& ReOmega) const
91 {
92  return(scalar(1) - exp(-sqr(max(ReOmega - CtsCrit_, scalar(0)))/Ats_));
93 }
94 
95 
96 tmp<volScalarField> kkLOmega::fTaul
97 (
98  const volScalarField& lambdaEff,
99  const volScalarField& ktL,
100  const volScalarField& Omega
101 ) const
102 {
103  return
104  (
105  scalar(1)
106  - exp
107  (
108  -CtauL_*ktL
109  /(
110  sqr(lambdaEff*Omega)
112  (
113  "rootVSmall",
115  rootVSmall
116  )
117  )
118  )
119  );
120 }
121 
122 
123 tmp<volScalarField> kkLOmega::alphaT
124 (
125  const volScalarField& lambdaEff,
126  const volScalarField& fv,
127  const volScalarField& ktS
128 ) const
129 {
130  return(fv*CmuStd_*sqrt(ktS)*lambdaEff);
131 }
132 
133 
134 tmp<volScalarField> kkLOmega::fOmega
135 (
136  const volScalarField& lambdaEff,
137  const volScalarField& lambdaT
138 ) const
139 {
140  return
141  (
142  scalar(1)
143  - exp
144  (
145  -0.41
146  *pow4
147  (
148  lambdaEff
149  / (
150  lambdaT
152  (
153  "ROTvSmall",
154  lambdaT.dimensions(),
155  rootVSmall
156  )
157  )
158  )
159  )
160  );
161 }
162 
163 
164 tmp<volScalarField> kkLOmega::phiBP(const volScalarField& Omega) const
165 {
166  return
167  (
168  min
169  (
170  max
171  (
172  kt_/nu()
173  / (
174  Omega
176  (
177  "ROTvSmall",
178  Omega.dimensions(),
179  rootVSmall
180  )
181  )
182  - CbpCrit_,
183  scalar(0)
184  ),
185  scalar(50)
186  )
187  );
188 }
189 
190 
191 tmp<volScalarField> kkLOmega::phiNAT
192 (
193  const volScalarField& ReOmega,
194  const volScalarField& fNatCrit
195 ) const
196 {
197  return
198  (
199  max
200  (
201  ReOmega
202  - CnatCrit_
203  / (
204  fNatCrit + dimensionedScalar(dimless, rootVSmall)
205  ),
206  scalar(0)
207  )
208  );
209 }
210 
211 
212 tmp<volScalarField> kkLOmega::D(const volScalarField& k) const
213 {
214  return nu()*magSqr(fvc::grad(sqrt(k)));
215 }
216 
217 
218 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
219 
221 {
222  // Currently this function is not implemented due to the complexity of
223  // evaluating nut. Better calculate nut at the end of correct()
225 }
226 
227 
228 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
229 
231 (
232  const geometricOneField& alpha,
233  const geometricOneField& rho,
234  const volVectorField& U,
235  const surfaceScalarField& alphaRhoPhi,
236  const surfaceScalarField& phi,
237  const viscosity& viscosity,
238  const word& type
239 )
240 :
241  eddyViscosity<incompressible::RASModel>
242  (
243  type,
244  alpha,
245  rho,
246  U,
247  alphaRhoPhi,
248  phi,
249  viscosity
250  ),
251 
252  A0_("A0", coeffDict(), 4.04),
253  As_("As", coeffDict(), 2.12),
254  Av_("Av", coeffDict(), 6.75),
255  Abp_("Abp", coeffDict(), 0.6),
256  Anat_("Anat", coeffDict(), 200),
257  Ats_("Ats", coeffDict(), 200),
258  CbpCrit_("CbpCrit", coeffDict(), 1.2),
259  Cnc_("Cnc", coeffDict(), 0.1),
260  CnatCrit_("CnatCrit", coeffDict(), 1250),
261  Cint_("Cint", coeffDict(), 0.75),
262  CtsCrit_("CtsCrit", coeffDict(), 1000),
263  CrNat_("CrNat", coeffDict(), 0.02),
264  C11_("C11", coeffDict(), 3.4e-6),
265  C12_("C12", coeffDict(), 1.0e-10),
266  CR_("CR", coeffDict(), 0.12),
267  CalphaTheta_("CalphaTheta", coeffDict(), 0.035),
268  Css_("Css", coeffDict(), 1.5),
269  CtauL_("CtauL", coeffDict(), 4360),
270  Cw1_("Cw1", coeffDict(), 0.44),
271  Cw2_("Cw2", coeffDict(), 0.92),
272  Cw3_("Cw3", coeffDict(), 0.3),
273  CwR_("CwR", coeffDict(), 1.5),
274  Clambda_("Clambda", coeffDict(), 2.495),
275  CmuStd_("CmuStd", coeffDict(), 0.09),
276  Prtheta_("Prtheta", coeffDict(), 0.85),
277  Sigmak_("Sigmak", coeffDict(), 1),
278  Sigmaw_("Sigmaw", coeffDict(), 1.17),
279  omegaMin_("omegaMin", dimless/dimTime, coeffDict(), small),
280  kt_
281  (
282  IOobject
283  (
284  this->groupName("kt"),
285  runTime_.name(),
286  mesh_,
289  ),
290  mesh_
291  ),
292  kl_
293  (
294  IOobject
295  (
296  this->groupName("kl"),
297  runTime_.name(),
298  mesh_,
301  ),
302  mesh_
303  ),
304  omega_
305  (
306  IOobject
307  (
308  this->groupName("omega"),
309  runTime_.name(),
310  mesh_,
313  ),
314  mesh_
315  ),
316  epsilon_
317  (
318  IOobject
319  (
320  "epsilon",
321  runTime_.name(),
322  mesh_
323  ),
324  kt_*omega_
325  )
326 {
327  bound(kt_, kMin_);
328  bound(kl_, kMin_);
329  bound(omega_, omegaMin_);
330  epsilon_ = kt_*omega_ + D(kl_) + D(kt_);
331 
332  if (type == typeName)
333  {
334  // Evaluating nut_ is complex so start from the field read from file
335  nut_.correctBoundaryConditions();
336  }
337 }
338 
339 
340 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
341 
342 bool kkLOmega::read()
343 {
345  {
346  A0_.readIfPresent(coeffDict());
347  As_.readIfPresent(coeffDict());
348  Av_.readIfPresent(coeffDict());
349  Abp_.readIfPresent(coeffDict());
350  Anat_.readIfPresent(coeffDict());
351  Abp_.readIfPresent(coeffDict());
352  Ats_.readIfPresent(coeffDict());
353  CbpCrit_.readIfPresent(coeffDict());
354  Cnc_.readIfPresent(coeffDict());
355  CnatCrit_.readIfPresent(coeffDict());
356  Cint_.readIfPresent(coeffDict());
357  CtsCrit_.readIfPresent(coeffDict());
358  CrNat_.readIfPresent(coeffDict());
359  C11_.readIfPresent(coeffDict());
360  C12_.readIfPresent(coeffDict());
361  CR_.readIfPresent(coeffDict());
362  CalphaTheta_.readIfPresent(coeffDict());
363  Css_.readIfPresent(coeffDict());
364  CtauL_.readIfPresent(coeffDict());
365  Cw1_.readIfPresent(coeffDict());
366  Cw2_.readIfPresent(coeffDict());
367  Cw3_.readIfPresent(coeffDict());
368  CwR_.readIfPresent(coeffDict());
369  Clambda_.readIfPresent(coeffDict());
370  CmuStd_.readIfPresent(coeffDict());
371  Prtheta_.readIfPresent(coeffDict());
372  Sigmak_.readIfPresent(coeffDict());
373  Sigmaw_.readIfPresent(coeffDict());
374 
375  return true;
376  }
377  else
378  {
379  return false;
380  }
381 }
382 
383 
384 void kkLOmega::validate()
385 {}
386 
387 
388 void kkLOmega::correct()
389 {
391 
392  if (!turbulence_)
393  {
394  return;
395  }
396 
397  volScalarField y(this->y());
398 
399  const volScalarField lambdaT(sqrt(kt_)/(omega_ + omegaMin_));
400 
401  const volScalarField lambdaEff(min(Clambda_*y, lambdaT));
402 
403  const volScalarField fw
404  (
405  pow
406  (
407  lambdaEff
408  /(lambdaT + dimensionedScalar(dimLength, rootVSmall)),
409  2.0/3.0
410  )
411  );
412 
413  tmp<volTensorField> tgradU(fvc::grad(U_));
414  const volTensorField& gradU = tgradU();
415 
416  const volScalarField Omega(sqrt(2.0)*mag(skew(gradU)));
417 
418  const volScalarField S2(2*magSqr(dev(symm(gradU))));
419 
420  const volScalarField ktS(fSS(Omega)*fw*kt_);
421 
422  const volScalarField nuts
423  (
424  fv(sqr(fw)*kt_/nu()/(omega_ + omegaMin_))
425  *fINT()
426  *Cmu(sqrt(S2))*sqrt(ktS)*lambdaEff
427  );
428  const volScalarField Pkt(this->GName(), nuts*S2);
429 
430  const volScalarField ktL(kt_ - ktS);
431  const volScalarField ReOmega(sqr(y)*Omega/nu());
432 
433  const volScalarField nutl
434  (
435  min
436  (
437  C11_*fTaul(lambdaEff, ktL, Omega)*Omega*sqr(lambdaEff)
438  *sqrt(ktL)*lambdaEff/nu()
439  + C12_*BetaTS(ReOmega)*sqr(sqr(lambdaEff/Clambda_)*Omega)/nu()
440  , 0.5*(kl_ + ktL)/(sqrt(S2) + omegaMin_)
441  )
442  );
443 
444  const volScalarField Pkl(nutl*S2);
445 
446  const volScalarField alphaTEff
447  (
448  alphaT(lambdaEff, fv(sqr(fw)*kt_/nu()/(omega_ + omegaMin_)), ktS)
449  );
450 
451  // By pass source term divided by kl_
452 
453  const dimensionedScalar fwMin("small", dimless, rootVSmall);
454 
455  const volScalarField Rbp
456  (
457  CR_*(1 - exp(-phiBP(Omega)()/Abp_))*omega_
458  /(fw + fwMin)
459  );
460 
461  const volScalarField fNatCrit(1 - exp(-Cnc_*sqrt(kl_)*y/nu()));
462 
463  // Natural source term divided by kl_
464  const volScalarField Rnat
465  (
466  CrNat_*(1 - exp(-phiNAT(ReOmega, fNatCrit)/Anat_))*Omega
467  );
468 
469 
471 
472  // Turbulence specific dissipation rate equation
473  tmp<fvScalarMatrix> omegaEqn
474  (
476  + fvm::div(phi_, omega_)
477  - fvm::laplacian(DomegaEff(alphaTEff), omega_)
478  ==
479  Cw1_*Pkt*omega_/(kt_ + kMin_)
480  - fvm::SuSp
481  (
482  (1 - CwR_/(fw + fwMin))*kl_*(Rbp + Rnat)/(kt_ + kMin_)
483  , omega_
484  )
485  - fvm::Sp(Cw2_*sqr(fw)*omega_, omega_)
486  + (
487  Cw3_*fOmega(lambdaEff, lambdaT)*alphaTEff*sqr(fw)*sqrt(kt_)
488  )()()/pow3(y())
489  );
490 
491  omegaEqn.ref().relax();
492  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
493 
494  solve(omegaEqn);
496 
497 
498  const volScalarField Dl(D(kl_));
499 
500  // Laminar kinetic energy equation
501  tmp<fvScalarMatrix> klEqn
502  (
503  fvm::ddt(kl_)
504  + fvm::div(phi_, kl_)
505  - fvm::laplacian(nu(), kl_)
506  ==
507  Pkl
508  - fvm::Sp(Rbp + Rnat + Dl/(kl_ + kMin_), kl_)
509  );
510 
511  klEqn.ref().relax();
512  klEqn.ref().boundaryManipulate(kl_.boundaryFieldRef());
513 
514  solve(klEqn);
515  bound(kl_, kMin_);
516 
517 
518  const volScalarField Dt(D(kt_));
519 
520  // Turbulent kinetic energy equation
521  tmp<fvScalarMatrix> ktEqn
522  (
523  fvm::ddt(kt_)
524  + fvm::div(phi_, kt_)
525  - fvm::laplacian(DkEff(alphaTEff), kt_)
526  ==
527  Pkt
528  + (Rbp + Rnat)*kl_
529  - fvm::Sp(omega_ + Dt/(kt_+ kMin_), kt_)
530  );
531 
532  ktEqn.ref().relax();
533  ktEqn.ref().boundaryManipulate(kt_.boundaryFieldRef());
534 
535  solve(ktEqn);
536  bound(kt_, kMin_);
537 
538 
539  // Update total fluctuation kinetic energy dissipation rate
540  epsilon_ = kt_*omega_ + Dl + Dt;
541 
542 
543  // Re-calculate turbulent viscosity
544  nut_ = nuts + nutl;
546 }
547 
548 
549 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
550 
551 } // End namespace RASModels
552 } // End namespace incompressible
553 } // End namespace Foam
554 
555 // ************************************************************************* //
scalar y
label k
Bound the given scalar field where it is below the specified minimum.
void updateCoeffs()
Update the boundary condition coefficients.
Boundary & boundaryFieldRef()
Return a reference to the boundary field.
void correctBoundaryConditions()
Correct boundary field.
bool readIfPresent(const dictionary &, const unitConversion &defaultUnits=NullObjectRef< unitConversion >())
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.
virtual void validate()
Validate the turbulence fields after construction.
virtual tmp< volScalarField > k() const
Return the total fluctuation kinetic energy.
Definition: kkLOmega.H:281
dimensionedScalar omegaMin_
Lower limit for omega.
Definition: kkLOmega.H:195
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
kkLOmega(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 bool read()
Read RASProperties dictionary.
tmp< volScalarField > DomegaEff(const volScalarField &alphaT) const
Return the effective diffusivity for omega.
Definition: kkLOmega.H:253
tmp< volScalarField > DkEff(const volScalarField &alphaT) const
Return the effective diffusivity for k.
Definition: kkLOmega.H:243
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:381
makeMomentumTransportModelTypes(geometricOneField, geometricOneField, incompressibleMomentumTransportModel)
Definition: kkLOmega.C:32
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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 > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
static const coefficient D("D", dimTemperature, 257.14)
Namespace for OpenFOAM.
void skew(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
dimensionedScalar exp(const dimensionedScalar &ds)
VolField< vector > volVectorField
Definition: volFieldsFwd.H:63
addToRunTimeSelectionTable(polyPatch, mergedCyclicPolyPatch, word)
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimLength
void pow4(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:66
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
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
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
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.
labelList fv(nPoints)
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))