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-2018 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 "fvOptions.H"
28 #include "bound.H"
29 #include "wallDist.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 
36 // * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * * //
37 
38 template<class TurbulenceModel, class BasicTurbulenceModel>
39 tmp<volScalarField>
41 (
42  const volScalarField& CDkOmega
43 ) const
44 {
45  tmp<volScalarField> CDkOmegaPlus = max
46  (
47  CDkOmega,
48  dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
49  );
50 
51  tmp<volScalarField> arg1 = min
52  (
53  min
54  (
55  max
56  (
57  (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
58  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
59  ),
60  (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
61  ),
62  scalar(10)
63  );
64 
65  return tanh(pow4(arg1));
66 }
67 
68 template<class TurbulenceModel, class BasicTurbulenceModel>
69 tmp<volScalarField>
71 {
72  tmp<volScalarField> arg2 = min
73  (
74  max
75  (
76  (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
77  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
78  ),
79  scalar(100)
80  );
81 
82  return tanh(sqr(arg2));
83 }
84 
85 template<class TurbulenceModel, class BasicTurbulenceModel>
86 tmp<volScalarField>
88 {
89  tmp<volScalarField> arg3 = min
90  (
91  150*(this->mu()/this->rho_)/(omega_*sqr(y_)),
92  scalar(10)
93  );
94 
95  return 1 - tanh(pow4(arg3));
96 }
97 
98 template<class TurbulenceModel, class BasicTurbulenceModel>
99 tmp<volScalarField>
100 kOmegaSST<TurbulenceModel, BasicTurbulenceModel>::kOmegaSST::F23() const
101 {
102  tmp<volScalarField> f23(F2());
103 
104  if (F3_)
105  {
106  f23.ref() *= F3();
107  }
108 
109  return f23;
110 }
111 
112 
113 template<class TurbulenceModel, class BasicTurbulenceModel>
115 (
116  const volScalarField& S2,
117  const volScalarField& F2
118 )
119 {
120  this->nut_ = a1_*k_/max(a1_*omega_, b1_*F2*sqrt(S2));
121  this->nut_.correctBoundaryConditions();
122  fv::options::New(this->mesh_).correct(this->nut_);
123 
124  BasicTurbulenceModel::correctNut();
125 }
126 
127 
128 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
129 
130 template<class TurbulenceModel, class BasicTurbulenceModel>
132 {
133  correctNut(2*magSqr(symm(fvc::grad(this->U_))), F23());
134 }
135 
136 
137 template<class TurbulenceModel, class BasicTurbulenceModel>
140 (
141  const volScalarField::Internal& G
142 ) const
143 {
144  return min(G, (c1_*betaStar_)*this->k_()*this->omega_());
145 }
146 
147 
148 template<class TurbulenceModel, class BasicTurbulenceModel>
151 (
152  const volScalarField::Internal& F1,
153  const volScalarField::Internal& F2
154 ) const
155 {
156  return betaStar_*omega_();
157 }
158 
159 
160 template<class TurbulenceModel, class BasicTurbulenceModel>
163 {
164  return tmp<fvScalarMatrix>
165  (
166  new fvScalarMatrix
167  (
168  k_,
169  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
170  )
171  );
172 }
173 
174 
175 template<class TurbulenceModel, class BasicTurbulenceModel>
178 {
179  return tmp<fvScalarMatrix>
180  (
181  new fvScalarMatrix
182  (
183  omega_,
184  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
185  )
186  );
187 }
188 
189 
190 template<class TurbulenceModel, class BasicTurbulenceModel>
192 (
193  const volScalarField::Internal& S2,
194  const volScalarField::Internal& gamma,
195  const volScalarField::Internal& beta
196 ) const
197 {
198  return tmp<fvScalarMatrix>
199  (
200  new fvScalarMatrix
201  (
202  omega_,
203  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
204  )
205  );
206 }
207 
208 
209 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
210 
211 template<class TurbulenceModel, class BasicTurbulenceModel>
213 (
214  const word& type,
215  const alphaField& alpha,
216  const rhoField& rho,
217  const volVectorField& U,
218  const surfaceScalarField& alphaRhoPhi,
219  const surfaceScalarField& phi,
220  const transportModel& transport,
221  const word& propertiesName
222 )
223 :
225  (
226  type,
227  alpha,
228  rho,
229  U,
230  alphaRhoPhi,
231  phi,
232  transport,
233  propertiesName
234  ),
235 
236  alphaK1_
237  (
239  (
240  "alphaK1",
241  this->coeffDict_,
242  0.85
243  )
244  ),
245  alphaK2_
246  (
248  (
249  "alphaK2",
250  this->coeffDict_,
251  1.0
252  )
253  ),
254  alphaOmega1_
255  (
257  (
258  "alphaOmega1",
259  this->coeffDict_,
260  0.5
261  )
262  ),
263  alphaOmega2_
264  (
266  (
267  "alphaOmega2",
268  this->coeffDict_,
269  0.856
270  )
271  ),
272  gamma1_
273  (
275  (
276  "gamma1",
277  this->coeffDict_,
278  5.0/9.0
279  )
280  ),
281  gamma2_
282  (
284  (
285  "gamma2",
286  this->coeffDict_,
287  0.44
288  )
289  ),
290  beta1_
291  (
293  (
294  "beta1",
295  this->coeffDict_,
296  0.075
297  )
298  ),
299  beta2_
300  (
302  (
303  "beta2",
304  this->coeffDict_,
305  0.0828
306  )
307  ),
308  betaStar_
309  (
311  (
312  "betaStar",
313  this->coeffDict_,
314  0.09
315  )
316  ),
317  a1_
318  (
320  (
321  "a1",
322  this->coeffDict_,
323  0.31
324  )
325  ),
326  b1_
327  (
329  (
330  "b1",
331  this->coeffDict_,
332  1.0
333  )
334  ),
335  c1_
336  (
338  (
339  "c1",
340  this->coeffDict_,
341  10.0
342  )
343  ),
344  F3_
345  (
347  (
348  "F3",
349  this->coeffDict_,
350  false
351  )
352  ),
353 
354  y_(wallDist::New(this->mesh_).y()),
355 
356  k_
357  (
358  IOobject
359  (
360  IOobject::groupName("k", alphaRhoPhi.group()),
361  this->runTime_.timeName(),
362  this->mesh_,
365  ),
366  this->mesh_
367  ),
368  omega_
369  (
370  IOobject
371  (
372  IOobject::groupName("omega", alphaRhoPhi.group()),
373  this->runTime_.timeName(),
374  this->mesh_,
377  ),
378  this->mesh_
379  )
380 {
381  bound(k_, this->kMin_);
382  bound(omega_, this->omegaMin_);
383 }
384 
385 
386 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
387 
388 template<class TurbulenceModel, class BasicTurbulenceModel>
390 {
391  if (TurbulenceModel::read())
392  {
393  alphaK1_.readIfPresent(this->coeffDict());
394  alphaK2_.readIfPresent(this->coeffDict());
395  alphaOmega1_.readIfPresent(this->coeffDict());
396  alphaOmega2_.readIfPresent(this->coeffDict());
397  gamma1_.readIfPresent(this->coeffDict());
398  gamma2_.readIfPresent(this->coeffDict());
399  beta1_.readIfPresent(this->coeffDict());
400  beta2_.readIfPresent(this->coeffDict());
401  betaStar_.readIfPresent(this->coeffDict());
402  a1_.readIfPresent(this->coeffDict());
403  b1_.readIfPresent(this->coeffDict());
404  c1_.readIfPresent(this->coeffDict());
405  F3_.readIfPresent("F3", this->coeffDict());
406 
407  return true;
408  }
409  else
410  {
411  return false;
412  }
413 }
414 
415 
416 template<class TurbulenceModel, class BasicTurbulenceModel>
418 {
419  if (!this->turbulence_)
420  {
421  return;
422  }
423 
424  // Local references
425  const alphaField& alpha = this->alpha_;
426  const rhoField& rho = this->rho_;
427  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
428  const volVectorField& U = this->U_;
429  volScalarField& nut = this->nut_;
430  fv::options& fvOptions(fv::options::New(this->mesh_));
431 
433 
435  (
436  fvc::div(fvc::absolute(this->phi(), U))()()
437  );
438 
439  tmp<volTensorField> tgradU = fvc::grad(U);
440  volScalarField S2(2*magSqr(symm(tgradU())));
441  volScalarField::Internal GbyNu(dev(twoSymm(tgradU()())) && tgradU()());
442  volScalarField::Internal G(this->GName(), nut()*GbyNu);
443  tgradU.clear();
444 
445  // Update omega and G at the wall
446  omega_.boundaryFieldRef().updateCoeffs();
447 
448  volScalarField CDkOmega
449  (
450  (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
451  );
452 
453  volScalarField F1(this->F1(CDkOmega));
454  volScalarField F23(this->F23());
455 
456  {
457  volScalarField::Internal gamma(this->gamma(F1));
459 
460  // Turbulent frequency equation
461  tmp<fvScalarMatrix> omegaEqn
462  (
463  fvm::ddt(alpha, rho, omega_)
464  + fvm::div(alphaRhoPhi, omega_)
465  - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
466  ==
467  alpha()*rho()*gamma
468  *min
469  (
470  GbyNu,
471  (c1_/a1_)*betaStar_*omega_()
472  *max(a1_*omega_(), b1_*F23()*sqrt(S2()))
473  )
474  - fvm::SuSp((2.0/3.0)*alpha()*rho()*gamma*divU, omega_)
475  - fvm::Sp(alpha()*rho()*beta*omega_(), omega_)
476  - fvm::SuSp
477  (
478  alpha()*rho()*(F1() - scalar(1))*CDkOmega()/omega_(),
479  omega_
480  )
481  + Qsas(S2(), gamma, beta)
482  + omegaSource()
483  + fvOptions(alpha, rho, omega_)
484  );
485 
486  omegaEqn.ref().relax();
487  fvOptions.constrain(omegaEqn.ref());
488  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
489  solve(omegaEqn);
490  fvOptions.correct(omega_);
491  bound(omega_, this->omegaMin_);
492  }
493 
494  // Turbulent kinetic energy equation
496  (
497  fvm::ddt(alpha, rho, k_)
498  + fvm::div(alphaRhoPhi, k_)
499  - fvm::laplacian(alpha*rho*DkEff(F1), k_)
500  ==
501  alpha()*rho()*Pk(G)
502  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
503  - fvm::Sp(alpha()*rho()*epsilonByk(F1, F23), k_)
504  + kSource()
505  + fvOptions(alpha, rho, k_)
506  );
507 
508  kEqn.ref().relax();
509  fvOptions.constrain(kEqn.ref());
510  solve(kEqn);
511  fvOptions.correct(k_);
512  bound(k_, this->kMin_);
513 
514  correctNut(S2, F23);
515 }
516 
517 
518 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
519 
520 } // End namespace Foam
521 
522 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:176
dimensionedScalar tanh(const dimensionedScalar &ds)
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
fv::options & fvOptions
surfaceScalarField & phi
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define F1(B, C, D)
Definition: SHA1.C:173
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
Templated abstract base class for turbulence models.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
dimensionedScalar sqrt(const dimensionedScalar &ds)
#define F2(B, C, D)
Definition: SHA1.C:174
virtual tmp< fvScalarMatrix > omegaSource() const
Generic dimensioned Type class.
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Finite-volume options.
Definition: fvOptions.H:52
static Switch lookupOrAddToDict(const word &, dictionary &, const Switch &defaultValue=false)
Construct from dictionary, supplying default value so that if the.
Definition: Switch.C:105
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
scalar y
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows...
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
#define F3(B, C, D)
Definition: SHA1.C:175
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
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
virtual tmp< fvScalarMatrix > Qsas(const volScalarField::Internal &S2, const volScalarField::Internal &gamma, const volScalarField::Internal &beta) const
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.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:72
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:523
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
const dimensionedScalar mu
Atomic mass unit.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
dimensionedScalar pow4(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< fvScalarMatrix > kSource() const
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:98
virtual tmp< volScalarField::Internal > epsilonByk(const volScalarField::Internal &F1, const volScalarField::Internal &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:691
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
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
virtual bool read()
Re-read model coefficients if they have changed.
scalar nut
Namespace for OpenFOAM.
virtual void correctNut()