kOmegaSSTBase.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 // * * * * * * * * * * * * Private 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>
139 (
140  const volScalarField& F1,
141  const volScalarField& F2
142 ) const
143 {
144  return betaStar_*omega_;
145 }
146 
147 
148 template<class TurbulenceModel, class BasicTurbulenceModel>
151 {
152  return tmp<fvScalarMatrix>
153  (
154  new fvScalarMatrix
155  (
156  k_,
157  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
158  )
159  );
160 }
161 
162 
163 template<class TurbulenceModel, class BasicTurbulenceModel>
166 {
167  return tmp<fvScalarMatrix>
168  (
169  new fvScalarMatrix
170  (
171  omega_,
172  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
173  )
174  );
175 }
176 
177 
178 template<class TurbulenceModel, class BasicTurbulenceModel>
180 (
181  const volScalarField& S2,
182  const volScalarField& gamma,
183  const volScalarField& beta
184 ) 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 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
198 
199 template<class TurbulenceModel, class BasicTurbulenceModel>
201 (
202  const word& type,
203  const alphaField& alpha,
204  const rhoField& rho,
205  const volVectorField& U,
206  const surfaceScalarField& alphaRhoPhi,
207  const surfaceScalarField& phi,
208  const transportModel& transport,
209  const word& propertiesName
210 )
211 :
213  (
214  type,
215  alpha,
216  rho,
217  U,
218  alphaRhoPhi,
219  phi,
220  transport,
221  propertiesName
222  ),
223 
224  alphaK1_
225  (
227  (
228  "alphaK1",
229  this->coeffDict_,
230  0.85
231  )
232  ),
233  alphaK2_
234  (
236  (
237  "alphaK2",
238  this->coeffDict_,
239  1.0
240  )
241  ),
242  alphaOmega1_
243  (
245  (
246  "alphaOmega1",
247  this->coeffDict_,
248  0.5
249  )
250  ),
251  alphaOmega2_
252  (
254  (
255  "alphaOmega2",
256  this->coeffDict_,
257  0.856
258  )
259  ),
260  gamma1_
261  (
263  (
264  "gamma1",
265  this->coeffDict_,
266  5.0/9.0
267  )
268  ),
269  gamma2_
270  (
272  (
273  "gamma2",
274  this->coeffDict_,
275  0.44
276  )
277  ),
278  beta1_
279  (
281  (
282  "beta1",
283  this->coeffDict_,
284  0.075
285  )
286  ),
287  beta2_
288  (
290  (
291  "beta2",
292  this->coeffDict_,
293  0.0828
294  )
295  ),
296  betaStar_
297  (
299  (
300  "betaStar",
301  this->coeffDict_,
302  0.09
303  )
304  ),
305  a1_
306  (
308  (
309  "a1",
310  this->coeffDict_,
311  0.31
312  )
313  ),
314  b1_
315  (
317  (
318  "b1",
319  this->coeffDict_,
320  1.0
321  )
322  ),
323  c1_
324  (
326  (
327  "c1",
328  this->coeffDict_,
329  10.0
330  )
331  ),
332  F3_
333  (
335  (
336  "F3",
337  this->coeffDict_,
338  false
339  )
340  ),
341 
342  y_(wallDist::New(this->mesh_).y()),
343 
344  k_
345  (
346  IOobject
347  (
348  IOobject::groupName("k", U.group()),
349  this->runTime_.timeName(),
350  this->mesh_,
353  ),
354  this->mesh_
355  ),
356  omega_
357  (
358  IOobject
359  (
360  IOobject::groupName("omega", U.group()),
361  this->runTime_.timeName(),
362  this->mesh_,
365  ),
366  this->mesh_
367  )
368 {
369  bound(k_, this->kMin_);
370  bound(omega_, this->omegaMin_);
371 }
372 
373 
374 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
375 
376 template<class TurbulenceModel, class BasicTurbulenceModel>
378 {
380  {
381  alphaK1_.readIfPresent(this->coeffDict());
382  alphaK2_.readIfPresent(this->coeffDict());
383  alphaOmega1_.readIfPresent(this->coeffDict());
384  alphaOmega2_.readIfPresent(this->coeffDict());
385  gamma1_.readIfPresent(this->coeffDict());
386  gamma2_.readIfPresent(this->coeffDict());
387  beta1_.readIfPresent(this->coeffDict());
388  beta2_.readIfPresent(this->coeffDict());
389  betaStar_.readIfPresent(this->coeffDict());
390  a1_.readIfPresent(this->coeffDict());
391  b1_.readIfPresent(this->coeffDict());
392  c1_.readIfPresent(this->coeffDict());
393  F3_.readIfPresent("F3", this->coeffDict());
394 
395  return true;
396  }
397  else
398  {
399  return false;
400  }
401 }
402 
403 
404 template<class TurbulenceModel, class BasicTurbulenceModel>
406 {
407  if (!this->turbulence_)
408  {
409  return;
410  }
411 
412  // Local references
413  const alphaField& alpha = this->alpha_;
414  const rhoField& rho = this->rho_;
415  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
416  const volVectorField& U = this->U_;
417  volScalarField& nut = this->nut_;
418  fv::options& fvOptions(fv::options::New(this->mesh_));
419 
421 
423 
424  tmp<volTensorField> tgradU = fvc::grad(U);
425  volScalarField S2(2*magSqr(symm(tgradU())));
426  volScalarField GbyNu((tgradU() && dev(twoSymm(tgradU()))));
427  volScalarField G(this->GName(), nut*GbyNu);
428  tgradU.clear();
429 
430  // Update omega and G at the wall
431  omega_.boundaryFieldRef().updateCoeffs();
432 
433  volScalarField CDkOmega
434  (
435  (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
436  );
437 
438  volScalarField F1(this->F1(CDkOmega));
439  volScalarField F23(this->F23());
440 
441  {
442  volScalarField gamma(this->gamma(F1));
443  volScalarField beta(this->beta(F1));
444 
445  // Turbulent frequency equation
446  tmp<fvScalarMatrix> omegaEqn
447  (
448  fvm::ddt(alpha, rho, omega_)
449  + fvm::div(alphaRhoPhi, omega_)
450  - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
451  ==
452  alpha*rho*gamma
453  *min
454  (
455  GbyNu,
456  (c1_/a1_)*betaStar_*omega_*max(a1_*omega_, b1_*F23*sqrt(S2))
457  )
458  - fvm::SuSp((2.0/3.0)*alpha*rho*gamma*divU, omega_)
459  - fvm::Sp(alpha*rho*beta*omega_, omega_)
460  - fvm::SuSp
461  (
462  alpha*rho*(F1 - scalar(1))*CDkOmega/omega_,
463  omega_
464  )
465  + Qsas(S2, gamma, beta)
466  + omegaSource()
467  + fvOptions(alpha, rho, omega_)
468  );
469 
470  omegaEqn.ref().relax();
471  fvOptions.constrain(omegaEqn.ref());
472  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
473  solve(omegaEqn);
474  fvOptions.correct(omega_);
475  bound(omega_, this->omegaMin_);
476  }
477 
478  // Turbulent kinetic energy equation
480  (
481  fvm::ddt(alpha, rho, k_)
482  + fvm::div(alphaRhoPhi, k_)
483  - fvm::laplacian(alpha*rho*DkEff(F1), k_)
484  ==
485  min(alpha*rho*G, (c1_*betaStar_)*alpha*rho*k_*omega_)
486  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
487  - fvm::Sp(alpha*rho*epsilonByk(F1, F23), k_)
488  + kSource()
489  + fvOptions(alpha, rho, k_)
490  );
491 
492  kEqn.ref().relax();
493  fvOptions.constrain(kEqn.ref());
494  solve(kEqn);
495  fvOptions.correct(k_);
496  bound(k_, this->kMin_);
497 
498  correctNut(S2, F23);
499 }
500 
501 
502 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
503 
504 } // End namespace Foam
505 
506 // ************************************************************************* //
surfaceScalarField & phi
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
const double e
Elementary charge.
Definition: doubleFloat.H:78
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
virtual tmp< fvScalarMatrix > kSource() const
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
BasicTurbulenceModel::alphaField alphaField
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
BasicTurbulenceModel::rhoField rhoField
virtual tmp< fvScalarMatrix > omegaSource() const
dimensionedScalar sqrt(const dimensionedScalar &ds)
#define F2(B, C, D)
Definition: SHA1.C:174
Generic dimensioned Type class.
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)
fv::options & fvOptions
static const wallDist & New(const fvMesh &mesh)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
scalar y
Implementation of the k-omega-SST turbulence model for incompressible and compressible flows...
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
#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)
virtual tmp< volScalarField > epsilonByk(const volScalarField &F1, const volScalarField &F2) const
Return epsilon/k which for standard RAS is betaStar*omega.
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 > &)
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-4.1/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:511
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
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
virtual tmp< fvScalarMatrix > Qsas(const volScalarField &S2, const volScalarField &gamma, const volScalarField &beta) const
dimensionedScalar pow4(const dimensionedScalar &ds)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
void boundaryManipulate(typename GeometricField< Type, fvPatchField, volMesh >::Boundary &values)
Manipulate based on a boundary field.
Definition: fvMatrix.C:679
A class for managing temporary objects.
Definition: PtrList.H:54
word group() const
Return group (extension part of name)
Definition: IOobject.C:239
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
virtual bool read()
Re-read model coefficients if they have changed.
scalar nut
Namespace for OpenFOAM.
BasicTurbulenceModel::transportModel transportModel
virtual void correctNut()