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-2015 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 "bound.H"
28 #include "wallDist.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 
35 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
36 
37 template<class BasicEddyViscosityModel>
39 (
40  const volScalarField& CDkOmega
41 ) const
42 {
43  tmp<volScalarField> CDkOmegaPlus = max
44  (
45  CDkOmega,
46  dimensionedScalar("1.0e-10", dimless/sqr(dimTime), 1.0e-10)
47  );
48 
50  (
51  min
52  (
53  max
54  (
55  (scalar(1)/betaStar_)*sqrt(k_)/(omega_*y_),
56  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
57  ),
58  (4*alphaOmega2_)*k_/(CDkOmegaPlus*sqr(y_))
59  ),
60  scalar(10)
61  );
62 
63  return tanh(pow4(arg1));
64 }
65 
66 
67 template<class BasicEddyViscosityModel>
69 {
71  (
72  max
73  (
74  (scalar(2)/betaStar_)*sqrt(k_)/(omega_*y_),
75  scalar(500)*(this->mu()/this->rho_)/(sqr(y_)*omega_)
76  ),
77  scalar(100)
78  );
79 
80  return tanh(sqr(arg2));
81 }
82 
83 
84 template<class BasicEddyViscosityModel>
86 {
88  (
89  150*(this->mu()/this->rho_)/(omega_*sqr(y_)),
90  scalar(10)
91  );
92 
93  return 1 - tanh(pow4(arg3));
94 }
95 
96 
97 template<class BasicEddyViscosityModel>
99 {
100  tmp<volScalarField> f23(F2());
101 
102  if (F3_)
103  {
104  f23() *= F3();
105  }
106 
107  return f23;
108 }
109 
110 
111 template<class BasicEddyViscosityModel>
113 (
114  const volScalarField& S2
115 )
116 {
117  // Correct the turbulence viscosity
118  this->nut_ = a1_*k_/max(a1_*omega_, b1_*F23()*sqrt(S2));
119  this->nut_.correctBoundaryConditions();
120 }
121 
122 
123 template<class BasicEddyViscosityModel>
125 {
126  correctNut(2*magSqr(symm(fvc::grad(this->U_))));
127 }
128 
129 
130 template<class BasicEddyViscosityModel>
132 {
133  return tmp<fvScalarMatrix>
134  (
135  new fvScalarMatrix
136  (
137  k_,
138  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
139  )
140  );
141 }
142 
143 
144 template<class BasicEddyViscosityModel>
146 {
147  return tmp<fvScalarMatrix>
148  (
149  new fvScalarMatrix
150  (
151  omega_,
152  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
153  )
154  );
155 }
156 
157 
158 template<class BasicEddyViscosityModel>
160 (
161  const volScalarField& S2,
162  const volScalarField& gamma,
163  const volScalarField& beta
164 ) const
165 {
166  return tmp<fvScalarMatrix>
167  (
168  new fvScalarMatrix
169  (
170  omega_,
171  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
172  )
173  );
174 }
175 
176 
177 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
178 
179 template<class BasicEddyViscosityModel>
181 (
182  const word& type,
183  const alphaField& alpha,
184  const rhoField& rho,
185  const volVectorField& U,
186  const surfaceScalarField& alphaRhoPhi,
187  const surfaceScalarField& phi,
188  const transportModel& transport,
189  const word& propertiesName
190 )
191 :
192  BasicEddyViscosityModel
193  (
194  type,
195  alpha,
196  rho,
197  U,
198  alphaRhoPhi,
199  phi,
200  transport,
201  propertiesName
202  ),
203 
204  alphaK1_
205  (
207  (
208  "alphaK1",
209  this->coeffDict_,
210  0.85
211  )
212  ),
213  alphaK2_
214  (
216  (
217  "alphaK2",
218  this->coeffDict_,
219  1.0
220  )
221  ),
222  alphaOmega1_
223  (
225  (
226  "alphaOmega1",
227  this->coeffDict_,
228  0.5
229  )
230  ),
231  alphaOmega2_
232  (
234  (
235  "alphaOmega2",
236  this->coeffDict_,
237  0.856
238  )
239  ),
240  gamma1_
241  (
243  (
244  "gamma1",
245  this->coeffDict_,
246  5.0/9.0
247  )
248  ),
249  gamma2_
250  (
252  (
253  "gamma2",
254  this->coeffDict_,
255  0.44
256  )
257  ),
258  beta1_
259  (
261  (
262  "beta1",
263  this->coeffDict_,
264  0.075
265  )
266  ),
267  beta2_
268  (
270  (
271  "beta2",
272  this->coeffDict_,
273  0.0828
274  )
275  ),
276  betaStar_
277  (
279  (
280  "betaStar",
281  this->coeffDict_,
282  0.09
283  )
284  ),
285  a1_
286  (
288  (
289  "a1",
290  this->coeffDict_,
291  0.31
292  )
293  ),
294  b1_
295  (
297  (
298  "b1",
299  this->coeffDict_,
300  1.0
301  )
302  ),
303  c1_
304  (
306  (
307  "c1",
308  this->coeffDict_,
309  10.0
310  )
311  ),
312  F3_
313  (
315  (
316  "F3",
317  this->coeffDict_,
318  false
319  )
320  ),
321 
322  y_(wallDist::New(this->mesh_).y()),
323 
324  k_
325  (
326  IOobject
327  (
328  IOobject::groupName("k", U.group()),
329  this->runTime_.timeName(),
330  this->mesh_,
333  ),
334  this->mesh_
335  ),
336  omega_
337  (
338  IOobject
339  (
340  IOobject::groupName("omega", U.group()),
341  this->runTime_.timeName(),
342  this->mesh_,
345  ),
346  this->mesh_
347  )
348 {
349  bound(k_, this->kMin_);
350  bound(omega_, this->omegaMin_);
351 }
352 
353 
354 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
355 
356 template<class BasicEddyViscosityModel>
358 {
360  {
361  alphaK1_.readIfPresent(this->coeffDict());
362  alphaK2_.readIfPresent(this->coeffDict());
363  alphaOmega1_.readIfPresent(this->coeffDict());
364  alphaOmega2_.readIfPresent(this->coeffDict());
365  gamma1_.readIfPresent(this->coeffDict());
366  gamma2_.readIfPresent(this->coeffDict());
367  beta1_.readIfPresent(this->coeffDict());
368  beta2_.readIfPresent(this->coeffDict());
369  betaStar_.readIfPresent(this->coeffDict());
370  a1_.readIfPresent(this->coeffDict());
371  b1_.readIfPresent(this->coeffDict());
372  c1_.readIfPresent(this->coeffDict());
373  F3_.readIfPresent("F3", this->coeffDict());
374 
375  return true;
376  }
377  else
378  {
379  return false;
380  }
381 }
382 
383 
384 template<class BasicEddyViscosityModel>
386 {
387  if (!this->turbulence_)
388  {
389  return;
390  }
391 
392  // Local references
393  const alphaField& alpha = this->alpha_;
394  const rhoField& rho = this->rho_;
395  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
396  const volVectorField& U = this->U_;
397  volScalarField& nut = this->nut_;
398 
400 
402 
403  tmp<volTensorField> tgradU = fvc::grad(U);
404  volScalarField S2(2*magSqr(symm(tgradU())));
405  volScalarField GbyNu((tgradU() && dev(twoSymm(tgradU()))));
406  volScalarField G(this->GName(), nut*GbyNu);
407  tgradU.clear();
408 
409  // Update omega and G at the wall
410  omega_.boundaryField().updateCoeffs();
411 
412  volScalarField CDkOmega
413  (
414  (2*alphaOmega2_)*(fvc::grad(k_) & fvc::grad(omega_))/omega_
415  );
416 
417  volScalarField F1(this->F1(CDkOmega));
418 
419  {
420  volScalarField gamma(this->gamma(F1));
421  volScalarField beta(this->beta(F1));
422 
423  // Turbulent frequency equation
424  tmp<fvScalarMatrix> omegaEqn
425  (
426  fvm::ddt(alpha, rho, omega_)
427  + fvm::div(alphaRhoPhi, omega_)
428  - fvm::laplacian(alpha*rho*DomegaEff(F1), omega_)
429  ==
430  alpha*rho*gamma
431  *min
432  (
433  GbyNu,
434  (c1_/a1_)*betaStar_*omega_*max(a1_*omega_, b1_*F23()*sqrt(S2))
435  )
436  - fvm::SuSp((2.0/3.0)*alpha*rho*gamma*divU, omega_)
437  - fvm::Sp(alpha*rho*beta*omega_, omega_)
438  - fvm::SuSp
439  (
440  alpha*rho*(F1 - scalar(1))*CDkOmega/omega_,
441  omega_
442  )
443  + Qsas(S2, gamma, beta)
444  + omegaSource()
445  );
446 
447  omegaEqn().relax();
448 
449  omegaEqn().boundaryManipulate(omega_.boundaryField());
450 
451  solve(omegaEqn);
452  bound(omega_, this->omegaMin_);
453  }
454 
455  // Turbulent kinetic energy equation
457  (
458  fvm::ddt(alpha, rho, k_)
459  + fvm::div(alphaRhoPhi, k_)
460  - fvm::laplacian(alpha*rho*DkEff(F1), k_)
461  ==
462  min(alpha*rho*G, (c1_*betaStar_)*alpha*rho*k_*omega_)
463  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
464  - fvm::Sp(alpha*rho*betaStar_*omega_, k_)
465  + kSource()
466  );
467 
468  kEqn().relax();
469  solve(kEqn);
470  bound(k_, this->kMin_);
471 
472  correctNut(S2);
473 }
474 
475 
476 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
477 
478 } // End namespace Foam
479 
480 // ************************************************************************* //
wallDist.H
Foam::kOmegaSSTBase::F2
tmp< volScalarField > F2() const
Definition: kOmegaSSTBase.C:68
Foam::kOmegaSSTBase::Qsas
virtual tmp< fvScalarMatrix > Qsas(const volScalarField &S2, const volScalarField &gamma, const volScalarField &beta) const
Definition: kOmegaSSTBase.C:160
beta
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
F1
#define F1(B, C, D)
Definition: SHA1.C:175
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::symm
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:82
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::constant::physicoChemical::mu
const dimensionedScalar mu
Atomic mass unit.
Definition: createFields.H:13
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::kOmegaSSTBase::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmegaSSTBase.C:131
Foam::kOmegaSSTBase::F1
tmp< volScalarField > F1(const volScalarField &CDkOmega) const
Definition: kOmegaSSTBase.C:39
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::kOmegaSSTBase::omegaSource
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmegaSSTBase.C:145
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::MeshObject< fvMesh, UpdateableMeshObject, wallDist >::New
static const wallDist & New(const fvMesh &mesh)
Definition: MeshObject.C:44
U
U
Definition: pEqn.H:46
Foam::kOmegaSSTBase::correctNut
virtual void correctNut()
Definition: kOmegaSSTBase.C:124
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:98
nut
nut
Definition: createTDFields.H:71
Foam::tanh
dimensionedScalar tanh(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:263
Foam::kOmegaSSTBase::kOmegaSSTBase
kOmegaSSTBase(const kOmegaSSTBase &)
kOmegaSSTBase.H
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
F3
#define F3(B, C, D)
Definition: SHA1.C:177
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::kOmegaSSTBase::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kOmegaSSTBase.C:357
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
Foam::e
const double e
Elementary charge.
Definition: doubleFloat.H:94
Foam::kOmegaSSTBase::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmegaSSTBase.C:385
rho
rho
Definition: pEqn.H:3
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::kOmegaSSTBase::F23
tmp< volScalarField > F23() const
Definition: kOmegaSSTBase.C:98
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
F2
#define F2(B, C, D)
Definition: SHA1.C:176
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
Foam::kOmegaSSTBase::F3
tmp< volScalarField > F3() const
Definition: kOmegaSSTBase.C:85
y
scalar y
Definition: LISASMDCalcMethod1.H:14
Foam::Switch::lookupOrAddToDict
static Switch lookupOrAddToDict(const word &, dictionary &, const Switch &defaultValue=false)
Construct from dictionary, supplying default value so that if the.
Definition: Switch.C:107
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104