kOmegaSSTBase.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 "kOmegaSSTBase.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "bound.H"
30 #include "wallDist.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 
37 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
38 
39 template<class MomentumTransportModel, class BasicMomentumTransportModel>
40 void
42 {
43  omega_ = max(omega_, k_/(this->nutMaxCoeff_*this->nu()));
44 }
45 
46 
47 template<class MomentumTransportModel, class BasicMomentumTransportModel>
50 (
51  const volScalarField& CDkOmega
52 ) const
53 {
54  tmp<volScalarField> CDkOmegaPlus = max
55  (
56  CDkOmega,
58  );
59 
61  (
62  min
63  (
64  max
65  (
66  (scalar(1)/betaStar_)*sqrt(k_)/(omega_*this->y()),
67  scalar(500)*this->nu()/(sqr(this->y())*omega_)
68  ),
69  (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(this->y()))
70  ),
71  scalar(10)
72  );
73 
74  return tanh(pow4(arg1));
75 }
76 
77 template<class MomentumTransportModel, class BasicMomentumTransportModel>
80 {
82  (
83  max
84  (
85  (scalar(2)/betaStar_)*sqrt(k_)/(omega_*this->y()),
86  scalar(500)*this->nu()/(sqr(this->y())*omega_)
87  ),
88  scalar(100)
89  );
90 
91  return tanh(sqr(arg2));
92 }
93 
94 template<class MomentumTransportModel, class BasicMomentumTransportModel>
97 {
99  (
100  150*this->nu()/(omega_*sqr(this->y())),
101  scalar(10)
102  );
103 
104  return 1 - tanh(pow4(arg3));
105 }
106 
107 template<class MomentumTransportModel, class BasicMomentumTransportModel>
110 {
111  tmp<volScalarField> f23(F2());
112 
113  if (F3_)
114  {
115  f23.ref() *= F3();
116  }
117 
118  return f23;
119 }
120 
121 
122 template<class MomentumTransportModel, class BasicMomentumTransportModel>
124 (
125  const volScalarField& S2,
126  const volScalarField& F2
127 )
128 {
129  this->nut_ = a1_*k_/max(a1_*omega_, b1_*F2*sqrt(S2));
130  this->nut_.correctBoundaryConditions();
131  fvConstraints::New(this->mesh_).constrain(this->nut_);
132 }
133 
134 
135 template<class MomentumTransportModel, class BasicMomentumTransportModel>
137 correctNut()
138 {
139  correctNut(2*magSqr(symm(fvc::grad(this->U_))), F23());
140 }
141 
142 
143 template<class MomentumTransportModel, class BasicMomentumTransportModel>
146 (
148 ) const
149 {
150  return min(G, (c1_*betaStar_)*this->k_()*this->omega_());
151 }
152 
153 
154 template<class MomentumTransportModel, class BasicMomentumTransportModel>
157 (
160 ) const
161 {
162  return betaStar_*omega_();
163 }
164 
165 
166 template<class MomentumTransportModel, class BasicMomentumTransportModel>
169 {
170  return tmp<fvScalarMatrix>
171  (
172  new fvScalarMatrix
173  (
174  k_,
175  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
176  )
177  );
178 }
179 
180 
181 template<class MomentumTransportModel, class BasicMomentumTransportModel>
184 omegaSource() const
185 {
186  return tmp<fvScalarMatrix>
187  (
188  new fvScalarMatrix
189  (
190  omega_,
191  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
192  )
193  );
194 }
195 
196 
197 template<class MomentumTransportModel, class BasicMomentumTransportModel>
200 (
201  const volScalarField::Internal& S2,
202  const volScalarField::Internal& gamma,
204 ) const
205 {
206  return tmp<fvScalarMatrix>
207  (
208  new fvScalarMatrix
209  (
210  omega_,
211  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
212  )
213  );
214 }
215 
216 
217 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
218 
219 template<class MomentumTransportModel, class BasicMomentumTransportModel>
221 (
222  const word& type,
223  const alphaField& alpha,
224  const rhoField& rho,
225  const volVectorField& U,
226  const surfaceScalarField& alphaRhoPhi,
227  const surfaceScalarField& phi,
229 )
230 :
231  MomentumTransportModel
232  (
233  type,
234  alpha,
235  rho,
236  U,
237  alphaRhoPhi,
238  phi,
239  viscosity
240  ),
241 
242  alphaK1_("alphaK1", this->coeffDict(), 0.85),
243  alphaK2_("alphaK2", this->coeffDict(), 1.0),
244  alphaOmega1_("alphaOmega1", this->coeffDict(), 0.5),
245  alphaOmega2_("alphaOmega2", this->coeffDict(), 0.856),
246  gamma1_("gamma1", this->coeffDict(), 5.0/9.0),
247  gamma2_("gamma2", this->coeffDict(), 0.44),
248  beta1_("beta1", this->coeffDict(), 0.075),
249  beta2_("beta2", this->coeffDict(), 0.0828),
250  betaStar_("betaStar", this->coeffDict(), 0.09),
251  a1_("a1", this->coeffDict(), 0.31),
252  b1_("b1", this->coeffDict(), 1.0),
253  c1_("c1", this->coeffDict(), 10.0),
254  F3_(this->coeffDict().template lookupOrDefault<Switch>("F3", false)),
255 
256  k_
257  (
258  IOobject
259  (
260  this->groupName("k"),
261  this->runTime_.name(),
262  this->mesh_,
263  IOobject::MUST_READ,
264  IOobject::AUTO_WRITE
265  ),
266  this->mesh_
267  ),
268  omega_
269  (
271  (
272  this->groupName("omega"),
273  this->runTime_.name(),
274  this->mesh_,
275  IOobject::MUST_READ,
276  IOobject::AUTO_WRITE
277  ),
278  this->mesh_
279  )
280 {
281  bound(k_, this->kMin_);
282  boundOmega();
283 }
284 
285 
286 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
287 
288 template<class MomentumTransportModel, class BasicMomentumTransportModel>
290 {
292  {
293  alphaK1_.readIfPresent(this->coeffDict());
294  alphaK2_.readIfPresent(this->coeffDict());
295  alphaOmega1_.readIfPresent(this->coeffDict());
296  alphaOmega2_.readIfPresent(this->coeffDict());
297  gamma1_.readIfPresent(this->coeffDict());
298  gamma2_.readIfPresent(this->coeffDict());
299  beta1_.readIfPresent(this->coeffDict());
300  beta2_.readIfPresent(this->coeffDict());
301  betaStar_.readIfPresent(this->coeffDict());
302  a1_.readIfPresent(this->coeffDict());
303  b1_.readIfPresent(this->coeffDict());
304  c1_.readIfPresent(this->coeffDict());
305  F3_.readIfPresent("F3", this->coeffDict());
306 
307  return true;
308  }
309  else
310  {
311  return false;
312  }
313 }
314 
315 
316 template<class MomentumTransportModel, class BasicMomentumTransportModel>
318 {
319  if (!this->turbulence_)
320  {
321  return;
322  }
323 
324  // Local references
325  const alphaField& alpha = this->alpha_;
326  const rhoField& rho = this->rho_;
327  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
328  const volVectorField& U = this->U_;
329  volScalarField& nut = this->nut_;
330  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
332  (
333  Foam::fvConstraints::New(this->mesh_)
334  );
335 
337 
339  (
340  fvc::div(fvc::absolute(this->phi(), U))()()
341  );
342 
343  tmp<volTensorField> tgradU = fvc::grad(U);
344  volScalarField S2(2*magSqr(symm(tgradU())));
345  volScalarField::Internal GbyNu(dev(twoSymm(tgradU()())) && tgradU()());
346  volScalarField::Internal G(this->GName(), nut()*GbyNu);
347  tgradU.clear();
348 
349  // Update omega and G at the wall
350  omega_.boundaryFieldRef().updateCoeffs();
351 
352  volScalarField CDkOmega
353  (
354  (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
355  );
356 
357  volScalarField F1(this->F1(CDkOmega));
358  volScalarField F23(this->F23());
359 
360  {
361  volScalarField::Internal gamma(this->gamma(F1));
362  volScalarField::Internal beta(this->beta(F1));
363 
364  // Turbulent frequency equation
365  tmp<fvScalarMatrix> omegaEqn
366  (
367  fvm::ddt(alpha, rho, omega_)
368  + fvm::div(alphaRhoPhi, omega_)
369  - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
370  ==
371  alpha()*rho()*gamma
372  *min
373  (
374  GbyNu,
375  (c1_/a1_)*betaStar_*omega_()
376  *max(a1_*omega_(), b1_*F23()*sqrt(S2()))
377  )
378  - fvm::SuSp((2.0/3.0)*alpha()*rho()*gamma*divU, omega_)
379  - fvm::Sp(alpha()*rho()*beta*omega_(), omega_)
380  - fvm::SuSp
381  (
382  alpha()*rho()*(F1() - scalar(1))*CDkOmega()/omega_(),
383  omega_
384  )
385  + Qsas(S2(), gamma, beta)
386  + omegaSource()
387  + fvModels.source(alpha, rho, omega_)
388  );
389 
390  omegaEqn.ref().relax();
391  fvConstraints.constrain(omegaEqn.ref());
392  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
393  solve(omegaEqn);
394  fvConstraints.constrain(omega_);
395  boundOmega();
396  }
397 
398  // Turbulent kinetic energy equation
400  (
401  fvm::ddt(alpha, rho, k_)
402  + fvm::div(alphaRhoPhi, k_)
403  - fvm::laplacian(alpha*rho*DkEff(F1), k_)
404  ==
405  alpha()*rho()*Pk(G)
406  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
407  - fvm::Sp(alpha()*rho()*epsilonByk(F1, F23), k_)
408  + kSource()
409  + fvModels.source(alpha, rho, k_)
410  );
411 
412  kEqn.ref().relax();
413  fvConstraints.constrain(kEqn.ref());
414  solve(kEqn);
416  bound(k_, this->kMin_);
417  boundOmega();
418 
419  correctNut(S2, F23);
420 }
421 
422 
423 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
424 
425 } // End namespace Foam
426 
427 // ************************************************************************* //
scalar y
#define F2(B, C, D)
Definition: SHA1.C:174
#define F1(B, C, D)
Definition: SHA1.C:173
#define F3(B, C, D)
Definition: SHA1.C:175
Bound the given scalar field where it is below the specified minimum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
void boundOmega()
Bound omega.
Definition: kOmegaSSTBase.C:41
virtual tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: kOmegaSSTBase.C:50
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
virtual tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:79
virtual tmp< volScalarField > F23() const
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
kOmegaSST(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity)
Construct from components.
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
virtual tmp< fvScalarMatrix > omegaSource() const
virtual void correctNut()
virtual tmp< fvScalarMatrix > kSource() const
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
virtual bool read()
Re-read model coefficients if they have changed.
virtual tmp< volScalarField > F3() const
Definition: kOmegaSSTBase.C:96
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:253
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
const scalar nut
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
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.
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimless
dimensionedScalar tanh(const dimensionedScalar &ds)
void pow4(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimVolume
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 sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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.