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-2020 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 MomentumTransportModel, class BasicMomentumTransportModel>
39 tmp<volScalarField>
41 (
42  const volScalarField& CDkOmega
43 ) const
44 {
45  tmp<volScalarField> CDkOmegaPlus = max
46  (
47  CDkOmega,
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 MomentumTransportModel, class BasicMomentumTransportModel>
69 tmp<volScalarField>
71 F2() const
72 {
73  tmp<volScalarField> arg2 = min
74  (
75  max
76  (
77  (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
78  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
79  ),
80  scalar(100)
81  );
82 
83  return tanh(sqr(arg2));
84 }
85 
86 template<class MomentumTransportModel, class BasicMomentumTransportModel>
87 tmp<volScalarField>
89 F3() const
90 {
91  tmp<volScalarField> arg3 = min
92  (
93  150*(this->mu()/this->rho_)/(omega_*sqr(y_)),
94  scalar(10)
95  );
96 
97  return 1 - tanh(pow4(arg3));
98 }
99 
100 template<class MomentumTransportModel, class BasicMomentumTransportModel>
101 tmp<volScalarField>
102 kOmegaSST<MomentumTransportModel, BasicMomentumTransportModel>::kOmegaSST::
103 F23() const
104 {
105  tmp<volScalarField> f23(F2());
106 
107  if (F3_)
108  {
109  f23.ref() *= F3();
110  }
111 
112  return f23;
113 }
114 
115 
116 template<class MomentumTransportModel, class BasicMomentumTransportModel>
118 (
119  const volScalarField& S2,
120  const volScalarField& F2
121 )
122 {
123  this->nut_ = a1_*k_/max(a1_*omega_, b1_*F2*sqrt(S2));
124  this->nut_.correctBoundaryConditions();
125  fv::options::New(this->mesh_).correct(this->nut_);
126 }
127 
128 
129 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
130 
131 template<class MomentumTransportModel, class BasicMomentumTransportModel>
134 {
135  correctNut(2*magSqr(symm(fvc::grad(this->U_))), F23());
136 }
137 
138 
139 template<class MomentumTransportModel, class BasicMomentumTransportModel>
142 (
143  const volScalarField::Internal& G
144 ) const
145 {
146  return min(G, (c1_*betaStar_)*this->k_()*this->omega_());
147 }
148 
149 
150 template<class MomentumTransportModel, class BasicMomentumTransportModel>
153 (
154  const volScalarField::Internal& F1,
155  const volScalarField::Internal& F2
156 ) const
157 {
158  return betaStar_*omega_();
159 }
160 
161 
162 template<class MomentumTransportModel, class BasicMomentumTransportModel>
165 {
166  return tmp<fvScalarMatrix>
167  (
168  new fvScalarMatrix
169  (
170  k_,
171  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
172  )
173  );
174 }
175 
176 
177 template<class MomentumTransportModel, class BasicMomentumTransportModel>
180 omegaSource() const
181 {
182  return tmp<fvScalarMatrix>
183  (
184  new fvScalarMatrix
185  (
186  omega_,
187  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
188  )
189  );
190 }
191 
192 
193 template<class MomentumTransportModel, class BasicMomentumTransportModel>
196 (
197  const volScalarField::Internal& S2,
198  const volScalarField::Internal& gamma,
199  const volScalarField::Internal& beta
200 ) const
201 {
202  return tmp<fvScalarMatrix>
203  (
204  new fvScalarMatrix
205  (
206  omega_,
207  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
208  )
209  );
210 }
211 
212 
213 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
214 
215 template<class MomentumTransportModel, class BasicMomentumTransportModel>
217 (
218  const word& type,
219  const alphaField& alpha,
220  const rhoField& rho,
221  const volVectorField& U,
222  const surfaceScalarField& alphaRhoPhi,
223  const surfaceScalarField& phi,
224  const transportModel& transport
225 )
226 :
228  (
229  type,
230  alpha,
231  rho,
232  U,
233  alphaRhoPhi,
234  phi,
235  transport
236  ),
237 
238  alphaK1_
239  (
241  (
242  "alphaK1",
243  this->coeffDict_,
244  0.85
245  )
246  ),
247  alphaK2_
248  (
250  (
251  "alphaK2",
252  this->coeffDict_,
253  1.0
254  )
255  ),
256  alphaOmega1_
257  (
259  (
260  "alphaOmega1",
261  this->coeffDict_,
262  0.5
263  )
264  ),
265  alphaOmega2_
266  (
268  (
269  "alphaOmega2",
270  this->coeffDict_,
271  0.856
272  )
273  ),
274  gamma1_
275  (
277  (
278  "gamma1",
279  this->coeffDict_,
280  5.0/9.0
281  )
282  ),
283  gamma2_
284  (
286  (
287  "gamma2",
288  this->coeffDict_,
289  0.44
290  )
291  ),
292  beta1_
293  (
295  (
296  "beta1",
297  this->coeffDict_,
298  0.075
299  )
300  ),
301  beta2_
302  (
304  (
305  "beta2",
306  this->coeffDict_,
307  0.0828
308  )
309  ),
310  betaStar_
311  (
313  (
314  "betaStar",
315  this->coeffDict_,
316  0.09
317  )
318  ),
319  a1_
320  (
322  (
323  "a1",
324  this->coeffDict_,
325  0.31
326  )
327  ),
328  b1_
329  (
331  (
332  "b1",
333  this->coeffDict_,
334  1.0
335  )
336  ),
337  c1_
338  (
340  (
341  "c1",
342  this->coeffDict_,
343  10.0
344  )
345  ),
346  F3_
347  (
349  (
350  "F3",
351  this->coeffDict_,
352  false
353  )
354  ),
355 
356  y_(wallDist::New(this->mesh_).y()),
357 
358  k_
359  (
360  IOobject
361  (
362  IOobject::groupName("k", alphaRhoPhi.group()),
363  this->runTime_.timeName(),
364  this->mesh_,
367  ),
368  this->mesh_
369  ),
370  omega_
371  (
372  IOobject
373  (
374  IOobject::groupName("omega", alphaRhoPhi.group()),
375  this->runTime_.timeName(),
376  this->mesh_,
379  ),
380  this->mesh_
381  )
382 {
383  bound(k_, this->kMin_);
384  bound(omega_, this->omegaMin_);
385 }
386 
387 
388 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
389 
390 template<class MomentumTransportModel, class BasicMomentumTransportModel>
392 {
394  {
395  alphaK1_.readIfPresent(this->coeffDict());
396  alphaK2_.readIfPresent(this->coeffDict());
397  alphaOmega1_.readIfPresent(this->coeffDict());
398  alphaOmega2_.readIfPresent(this->coeffDict());
399  gamma1_.readIfPresent(this->coeffDict());
400  gamma2_.readIfPresent(this->coeffDict());
401  beta1_.readIfPresent(this->coeffDict());
402  beta2_.readIfPresent(this->coeffDict());
403  betaStar_.readIfPresent(this->coeffDict());
404  a1_.readIfPresent(this->coeffDict());
405  b1_.readIfPresent(this->coeffDict());
406  c1_.readIfPresent(this->coeffDict());
407  F3_.readIfPresent("F3", this->coeffDict());
408 
409  return true;
410  }
411  else
412  {
413  return false;
414  }
415 }
416 
417 
418 template<class MomentumTransportModel, class BasicMomentumTransportModel>
420 {
421  if (!this->turbulence_)
422  {
423  return;
424  }
425 
426  // Local references
427  const alphaField& alpha = this->alpha_;
428  const rhoField& rho = this->rho_;
429  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
430  const volVectorField& U = this->U_;
431  volScalarField& nut = this->nut_;
432  fv::options& fvOptions(fv::options::New(this->mesh_));
433 
435 
437  (
438  fvc::div(fvc::absolute(this->phi(), U))()()
439  );
440 
441  tmp<volTensorField> tgradU = fvc::grad(U);
442  volScalarField S2(2*magSqr(symm(tgradU())));
443  volScalarField::Internal GbyNu(dev(twoSymm(tgradU()())) && tgradU()());
444  volScalarField::Internal G(this->GName(), nut()*GbyNu);
445  tgradU.clear();
446 
447  // Update omega and G at the wall
448  omega_.boundaryFieldRef().updateCoeffs();
449 
450  volScalarField CDkOmega
451  (
452  (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
453  );
454 
455  volScalarField F1(this->F1(CDkOmega));
456  volScalarField F23(this->F23());
457 
458  {
459  volScalarField::Internal gamma(this->gamma(F1));
460  volScalarField::Internal beta(this->beta(F1));
461 
462  // Turbulent frequency equation
463  tmp<fvScalarMatrix> omegaEqn
464  (
465  fvm::ddt(alpha, rho, omega_)
466  + fvm::div(alphaRhoPhi, omega_)
467  - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
468  ==
469  alpha()*rho()*gamma
470  *min
471  (
472  GbyNu,
473  (c1_/a1_)*betaStar_*omega_()
474  *max(a1_*omega_(), b1_*F23()*sqrt(S2()))
475  )
476  - fvm::SuSp((2.0/3.0)*alpha()*rho()*gamma*divU, omega_)
477  - fvm::Sp(alpha()*rho()*beta*omega_(), omega_)
478  - fvm::SuSp
479  (
480  alpha()*rho()*(F1() - scalar(1))*CDkOmega()/omega_(),
481  omega_
482  )
483  + Qsas(S2(), gamma, beta)
484  + omegaSource()
485  + fvOptions(alpha, rho, omega_)
486  );
487 
488  omegaEqn.ref().relax();
489  fvOptions.constrain(omegaEqn.ref());
490  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
491  solve(omegaEqn);
492  fvOptions.correct(omega_);
493  bound(omega_, this->omegaMin_);
494  }
495 
496  // Turbulent kinetic energy equation
498  (
499  fvm::ddt(alpha, rho, k_)
500  + fvm::div(alphaRhoPhi, k_)
501  - fvm::laplacian(alpha*rho*DkEff(F1), k_)
502  ==
503  alpha()*rho()*Pk(G)
504  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
505  - fvm::Sp(alpha()*rho()*epsilonByk(F1, F23), k_)
506  + kSource()
507  + fvOptions(alpha, rho, k_)
508  );
509 
510  kEqn.ref().relax();
511  fvOptions.constrain(kEqn.ref());
512  solve(kEqn);
513  fvOptions.correct(k_);
514  bound(k_, this->kMin_);
515 
516  correctNut(S2, F23);
517 }
518 
519 
520 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
521 
522 } // End namespace Foam
523 
524 // ************************************************************************* //
static word group(const word &name)
Return group (extension part of name)
Definition: IOobject.C:144
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
kOmegaSST(const word &type, const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport)
Construct from components.
fv::options & fvOptions
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.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
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
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:88
static Switch lookupOrAddToDict(const word &, dictionary &, const Switch &defaultValue=false)
Construct from dictionary, supplying default value so that if the.
Definition: Switch.C:111
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
scalar y
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
BasicMomentumTransportModel::transportModel transportModel
Definition: RASModel.H:90
#define F3(B, C, D)
Definition: SHA1.C:175
A class for handling words, derived from string.
Definition: word.H:59
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(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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:68
void relax(const scalar alpha)
Relax matrix (for steady-state solution).
Definition: fvMatrix.C:521
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
virtual tmp< volScalarField::Internal > Pk(const volScalarField::Internal &G) const
Return k production rate.
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:89
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
const dimensionedScalar & mu
Atomic mass unit.
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
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
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:692
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
Templated abstract base class for turbulence models.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
virtual bool read()
Re-read model coefficients if they have changed.
scalar nut
Namespace for OpenFOAM.
virtual void correctNut()