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-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 "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_
253  (
255  (
256  "A0",
257  coeffDict_,
258  4.04
259  )
260  ),
261  As_
262  (
264  (
265  "As",
266  coeffDict_,
267  2.12
268  )
269  ),
270  Av_
271  (
273  (
274  "Av",
275  coeffDict_,
276  6.75
277  )
278  ),
279  Abp_
280  (
282  (
283  "Abp",
284  coeffDict_,
285  0.6
286  )
287  ),
288  Anat_
289  (
291  (
292  "Anat",
293  coeffDict_,
294  200
295  )
296  ),
297  Ats_
298  (
300  (
301  "Ats",
302  coeffDict_,
303  200
304  )
305  ),
306  CbpCrit_
307  (
309  (
310  "CbpCrit",
311  coeffDict_,
312  1.2
313  )
314  ),
315  Cnc_
316  (
318  (
319  "Cnc",
320  coeffDict_,
321  0.1
322  )
323  ),
324  CnatCrit_
325  (
327  (
328  "CnatCrit",
329  coeffDict_,
330  1250
331  )
332  ),
333  Cint_
334  (
336  (
337  "Cint",
338  coeffDict_,
339  0.75
340  )
341  ),
342  CtsCrit_
343  (
345  (
346  "CtsCrit",
347  coeffDict_,
348  1000
349  )
350  ),
351  CrNat_
352  (
354  (
355  "CrNat",
356  coeffDict_,
357  0.02
358  )
359  ),
360  C11_
361  (
363  (
364  "C11",
365  coeffDict_,
366  3.4e-6
367  )
368  ),
369  C12_
370  (
372  (
373  "C12",
374  coeffDict_,
375  1.0e-10
376  )
377  ),
378  CR_
379  (
381  (
382  "CR",
383  coeffDict_,
384  0.12
385  )
386  ),
387  CalphaTheta_
388  (
390  (
391  "CalphaTheta",
392  coeffDict_,
393  0.035
394  )
395  ),
396  Css_
397  (
399  (
400  "Css",
401  coeffDict_,
402  1.5
403  )
404  ),
405  CtauL_
406  (
408  (
409  "CtauL",
410  coeffDict_,
411  4360
412  )
413  ),
414  Cw1_
415  (
417  (
418  "Cw1",
419  coeffDict_,
420  0.44
421  )
422  ),
423  Cw2_
424  (
426  (
427  "Cw2",
428  coeffDict_,
429  0.92
430  )
431  ),
432  Cw3_
433  (
435  (
436  "Cw3",
437  coeffDict_,
438  0.3
439  )
440  ),
441  CwR_
442  (
444  (
445  "CwR",
446  coeffDict_,
447  1.5
448  )
449  ),
450  Clambda_
451  (
453  (
454  "Clambda",
455  coeffDict_,
456  2.495
457  )
458  ),
459  CmuStd_
460  (
462  (
463  "CmuStd",
464  coeffDict_,
465  0.09
466  )
467  ),
468  Prtheta_
469  (
471  (
472  "Prtheta",
473  coeffDict_,
474  0.85
475  )
476  ),
477  Sigmak_
478  (
480  (
481  "Sigmak",
482  coeffDict_,
483  1
484  )
485  ),
486  Sigmaw_
487  (
489  (
490  "Sigmaw",
491  coeffDict_,
492  1.17
493  )
494  ),
495  omegaMin_
496  (
498  (
499  "omegaMin",
500  coeffDict_,
502  small
503  )
504  ),
505  kt_
506  (
507  IOobject
508  (
509  this->groupName("kt"),
510  runTime_.name(),
511  mesh_,
514  ),
515  mesh_
516  ),
517  kl_
518  (
519  IOobject
520  (
521  this->groupName("kl"),
522  runTime_.name(),
523  mesh_,
526  ),
527  mesh_
528  ),
529  omega_
530  (
531  IOobject
532  (
533  this->groupName("omega"),
534  runTime_.name(),
535  mesh_,
538  ),
539  mesh_
540  ),
541  epsilon_
542  (
543  IOobject
544  (
545  "epsilon",
546  runTime_.name(),
547  mesh_
548  ),
549  kt_*omega_
550  )
551 {
552  bound(kt_, kMin_);
553  bound(kl_, kMin_);
554  bound(omega_, omegaMin_);
555  epsilon_ = kt_*omega_ + D(kl_) + D(kt_);
556 
557  if (type == typeName)
558  {
559  // Evaluating nut_ is complex so start from the field read from file
560  nut_.correctBoundaryConditions();
561 
562  printCoeffs(type);
563  }
564 }
565 
566 
567 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
568 
569 bool kkLOmega::read()
570 {
572  {
573  A0_.readIfPresent(coeffDict());
574  As_.readIfPresent(coeffDict());
575  Av_.readIfPresent(coeffDict());
576  Abp_.readIfPresent(coeffDict());
577  Anat_.readIfPresent(coeffDict());
578  Abp_.readIfPresent(coeffDict());
579  Ats_.readIfPresent(coeffDict());
580  CbpCrit_.readIfPresent(coeffDict());
581  Cnc_.readIfPresent(coeffDict());
582  CnatCrit_.readIfPresent(coeffDict());
583  Cint_.readIfPresent(coeffDict());
584  CtsCrit_.readIfPresent(coeffDict());
585  CrNat_.readIfPresent(coeffDict());
586  C11_.readIfPresent(coeffDict());
587  C12_.readIfPresent(coeffDict());
588  CR_.readIfPresent(coeffDict());
589  CalphaTheta_.readIfPresent(coeffDict());
590  Css_.readIfPresent(coeffDict());
591  CtauL_.readIfPresent(coeffDict());
592  Cw1_.readIfPresent(coeffDict());
593  Cw2_.readIfPresent(coeffDict());
594  Cw3_.readIfPresent(coeffDict());
595  CwR_.readIfPresent(coeffDict());
596  Clambda_.readIfPresent(coeffDict());
597  CmuStd_.readIfPresent(coeffDict());
598  Prtheta_.readIfPresent(coeffDict());
599  Sigmak_.readIfPresent(coeffDict());
600  Sigmaw_.readIfPresent(coeffDict());
601 
602  return true;
603  }
604  else
605  {
606  return false;
607  }
608 }
609 
610 
611 void kkLOmega::validate()
612 {}
613 
614 
615 void kkLOmega::correct()
616 {
618 
619  if (!turbulence_)
620  {
621  return;
622  }
623 
624  volScalarField y(this->y());
625 
626  const volScalarField lambdaT(sqrt(kt_)/(omega_ + omegaMin_));
627 
628  const volScalarField lambdaEff(min(Clambda_*y, lambdaT));
629 
630  const volScalarField fw
631  (
632  pow
633  (
634  lambdaEff
635  /(lambdaT + dimensionedScalar(dimLength, rootVSmall)),
636  2.0/3.0
637  )
638  );
639 
640  tmp<volTensorField> tgradU(fvc::grad(U_));
641  const volTensorField& gradU = tgradU();
642 
643  const volScalarField Omega(sqrt(2.0)*mag(skew(gradU)));
644 
645  const volScalarField S2(2*magSqr(dev(symm(gradU))));
646 
647  const volScalarField ktS(fSS(Omega)*fw*kt_);
648 
649  const volScalarField nuts
650  (
651  fv(sqr(fw)*kt_/nu()/(omega_ + omegaMin_))
652  *fINT()
653  *Cmu(sqrt(S2))*sqrt(ktS)*lambdaEff
654  );
655  const volScalarField Pkt(this->GName(), nuts*S2);
656 
657  const volScalarField ktL(kt_ - ktS);
658  const volScalarField ReOmega(sqr(y)*Omega/nu());
659 
660  const volScalarField nutl
661  (
662  min
663  (
664  C11_*fTaul(lambdaEff, ktL, Omega)*Omega*sqr(lambdaEff)
665  *sqrt(ktL)*lambdaEff/nu()
666  + C12_*BetaTS(ReOmega)*sqr(sqr(lambdaEff/Clambda_)*Omega)/nu()
667  , 0.5*(kl_ + ktL)/(sqrt(S2) + omegaMin_)
668  )
669  );
670 
671  const volScalarField Pkl(nutl*S2);
672 
673  const volScalarField alphaTEff
674  (
675  alphaT(lambdaEff, fv(sqr(fw)*kt_/nu()/(omega_ + omegaMin_)), ktS)
676  );
677 
678  // By pass source term divided by kl_
679 
680  const dimensionedScalar fwMin("small", dimless, rootVSmall);
681 
682  const volScalarField Rbp
683  (
684  CR_*(1 - exp(-phiBP(Omega)()/Abp_))*omega_
685  /(fw + fwMin)
686  );
687 
688  const volScalarField fNatCrit(1 - exp(-Cnc_*sqrt(kl_)*y/nu()));
689 
690  // Natural source term divided by kl_
691  const volScalarField Rnat
692  (
693  CrNat_*(1 - exp(-phiNAT(ReOmega, fNatCrit)/Anat_))*Omega
694  );
695 
696 
698 
699  // Turbulence specific dissipation rate equation
700  tmp<fvScalarMatrix> omegaEqn
701  (
703  + fvm::div(phi_, omega_)
704  - fvm::laplacian(DomegaEff(alphaTEff), omega_)
705  ==
706  Cw1_*Pkt*omega_/(kt_ + kMin_)
707  - fvm::SuSp
708  (
709  (1 - CwR_/(fw + fwMin))*kl_*(Rbp + Rnat)/(kt_ + kMin_)
710  , omega_
711  )
712  - fvm::Sp(Cw2_*sqr(fw)*omega_, omega_)
713  + (
714  Cw3_*fOmega(lambdaEff, lambdaT)*alphaTEff*sqr(fw)*sqrt(kt_)
715  )()()/pow3(y())
716  );
717 
718  omegaEqn.ref().relax();
719  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
720 
721  solve(omegaEqn);
723 
724 
725  const volScalarField Dl(D(kl_));
726 
727  // Laminar kinetic energy equation
728  tmp<fvScalarMatrix> klEqn
729  (
730  fvm::ddt(kl_)
731  + fvm::div(phi_, kl_)
732  - fvm::laplacian(nu(), kl_)
733  ==
734  Pkl
735  - fvm::Sp(Rbp + Rnat + Dl/(kl_ + kMin_), kl_)
736  );
737 
738  klEqn.ref().relax();
739  klEqn.ref().boundaryManipulate(kl_.boundaryFieldRef());
740 
741  solve(klEqn);
742  bound(kl_, kMin_);
743 
744 
745  const volScalarField Dt(D(kt_));
746 
747  // Turbulent kinetic energy equation
748  tmp<fvScalarMatrix> ktEqn
749  (
750  fvm::ddt(kt_)
751  + fvm::div(phi_, kt_)
752  - fvm::laplacian(DkEff(alphaTEff), kt_)
753  ==
754  Pkt
755  + (Rbp + Rnat)*kl_
756  - fvm::Sp(omega_ + Dt/(kt_+ kMin_), kt_)
757  );
758 
759  ktEqn.ref().relax();
760  ktEqn.ref().boundaryManipulate(kt_.boundaryFieldRef());
761 
762  solve(ktEqn);
763  bound(kt_, kMin_);
764 
765 
766  // Update total fluctuation kinetic energy dissipation rate
767  epsilon_ = kt_*omega_ + Dl + Dt;
768 
769 
770  // Re-calculate turbulent viscosity
771  nut_ = nuts + nutl;
773 }
774 
775 
776 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
777 
778 } // End namespace RASModels
779 } // End namespace incompressible
780 } // End namespace Foam
781 
782 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
scalar y
label k
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.
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
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)
const dimensionSet dimless
SurfaceField< scalar > surfaceScalarField
const dimensionSet dimLength
VolField< tensor > volTensorField
Definition: volFieldsFwd.H:68
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:64
dimensioned< scalar > mag(const dimensioned< Type > &)
defineTypeNameAndDebug(combustionModel, 0)
dimensionedScalar pow4(const dimensionedScalar &ds)
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
dimensioned< scalar > magSqr(const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
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))