kOmega2006.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) 2021-2023 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 "kOmega2006.H"
27 #include "fvModels.H"
28 #include "fvConstraints.H"
29 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
42 {
43  omega_ = max(omega_, k_/(this->nutMaxCoeff_*this->nu()));
44 }
45 
46 
47 template<class BasicMomentumTransportModel>
49 (
50  const volTensorField& gradU
51 )
52 {
53  this->nut_ = k_/max(omega_, Clim_*sqrt(2/betaStar_)*mag(dev(symm(gradU))));
54  this->nut_.correctBoundaryConditions();
55  fvConstraints::New(this->mesh_).constrain(this->nut_);
56 }
57 
58 
59 template<class BasicMomentumTransportModel>
60 tmp<volScalarField::Internal> kOmega2006<BasicMomentumTransportModel>::beta
61 (
62  const volTensorField& gradU
63 ) const
64 {
65  const volSymmTensorField::Internal S(symm(gradU()));
66  const volSymmTensorField::Internal Shat(S - 0.5*tr(S)*I);
67  const volTensorField::Internal Omega(skew(gradU.v()));
68 
69  const volScalarField::Internal ChiOmega
70  (
71  typedName("ChiOmega"),
72  mag((Omega & Omega) && Shat)/pow3(betaStar_*omega_.v())
73  );
74 
75  const volScalarField::Internal fBeta
76  (
77  typedName("fBeta"),
78  (1 + 85*ChiOmega)/(1 + 100*ChiOmega)
79  );
80 
81  return beta0_*fBeta;
82 }
83 
84 
85 template<class BasicMomentumTransportModel>
86 tmp<volScalarField::Internal>
87 kOmega2006<BasicMomentumTransportModel>::CDkOmega() const
88 {
89  return max
90  (
91  sigmaDo_*(fvc::grad(k_)().v() & fvc::grad(omega_)().v())/omega_(),
93  );
94 }
95 
96 
97 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
98 
99 template<class BasicMomentumTransportModel>
101 {
102  correctNut(fvc::grad(this->U_));
103 }
104 
105 
106 template<class BasicMomentumTransportModel>
108 {
109  return tmp<fvScalarMatrix>
110  (
111  new fvScalarMatrix
112  (
113  k_,
114  dimVolume*this->rho_.dimensions()*k_.dimensions()
115  /dimTime
116  )
117  );
118 }
119 
120 
121 template<class BasicMomentumTransportModel>
123 {
124  return tmp<fvScalarMatrix>
125  (
126  new fvScalarMatrix
127  (
128  omega_,
129  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
130  )
131  );
132 }
133 
134 
135 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
136 
137 template<class BasicMomentumTransportModel>
139 (
140  const alphaField& alpha,
141  const rhoField& rho,
142  const volVectorField& U,
143  const surfaceScalarField& alphaRhoPhi,
144  const surfaceScalarField& phi,
145  const viscosity& viscosity,
146  const word& type
147 )
148 :
149  eddyViscosity<RASModel<BasicMomentumTransportModel>>
150  (
151  type,
152  alpha,
153  rho,
154  U,
155  alphaRhoPhi,
156  phi,
157  viscosity
158  ),
159 
160  betaStar_
161  (
162  dimensioned<scalar>::lookupOrAddToDict
163  (
164  "betaStar",
165  this->coeffDict_,
166  0.09
167  )
168  ),
169  beta0_
170  (
171  dimensioned<scalar>::lookupOrAddToDict
172  (
173  "beta0",
174  this->coeffDict_,
175  0.0708
176  )
177  ),
178  gamma_
179  (
180  dimensioned<scalar>::lookupOrAddToDict
181  (
182  "gamma",
183  this->coeffDict_,
184  0.52
185  )
186  ),
187  Clim_
188  (
189  dimensioned<scalar>::lookupOrAddToDict
190  (
191  "Clim",
192  this->coeffDict_,
193  0.875
194  )
195  ),
196  sigmaDo_
197  (
198  dimensioned<scalar>::lookupOrAddToDict
199  (
200  "sigmaDo",
201  this->coeffDict_,
202  0.125
203  )
204  ),
205  alphaK_
206  (
207  dimensioned<scalar>::lookupOrAddToDict
208  (
209  "alphaK",
210  this->coeffDict_,
211  0.6
212  )
213  ),
214  alphaOmega_
215  (
216  dimensioned<scalar>::lookupOrAddToDict
217  (
218  "alphaOmega",
219  this->coeffDict_,
220  0.5
221  )
222  ),
223 
224  k_
225  (
226  IOobject
227  (
228  this->groupName("k"),
229  this->runTime_.name(),
230  this->mesh_,
231  IOobject::MUST_READ,
232  IOobject::AUTO_WRITE
233  ),
234  this->mesh_
235  ),
236  omega_
237  (
238  IOobject
239  (
240  this->groupName("omega"),
241  this->runTime_.name(),
242  this->mesh_,
243  IOobject::MUST_READ,
244  IOobject::AUTO_WRITE
245  ),
246  this->mesh_
247  )
248 {
249  bound(k_, this->kMin_);
250  boundOmega();
251 
252  if (type == typeName)
253  {
254  this->printCoeffs(type);
255  }
256 }
257 
258 
259 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
260 
261 template<class BasicMomentumTransportModel>
263 {
265  {
266  betaStar_.readIfPresent(this->coeffDict());
267  beta0_.readIfPresent(this->coeffDict());
268  gamma_.readIfPresent(this->coeffDict());
269  sigmaDo_.readIfPresent(this->coeffDict());
270  alphaK_.readIfPresent(this->coeffDict());
271  alphaOmega_.readIfPresent(this->coeffDict());
272 
273  return true;
274  }
275  else
276  {
277  return false;
278  }
279 }
280 
281 
282 template<class BasicMomentumTransportModel>
284 {
285  if (!this->turbulence_)
286  {
287  return;
288  }
289 
290  // Local references
291  const alphaField& alpha = this->alpha_;
292  const rhoField& rho = this->rho_;
293  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
294  const volVectorField& U = this->U_;
295  volScalarField& nut = this->nut_;
296  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
298  (
299  Foam::fvConstraints::New(this->mesh_)
300  );
301 
303 
305  (
306  typedName("divU"),
307  fvc::div(fvc::absolute(this->phi(), U))().v()
308  );
309 
310  const volTensorField gradU(fvc::grad(U));
312  (
313  this->GName(),
314  nut.v()*(dev(twoSymm(gradU.v())) && gradU.v())
315  );
316 
317  // Update omega and G at the wall
318  omega_.boundaryFieldRef().updateCoeffs();
319 
320  // Turbulence specific dissipation rate equation
321  tmp<fvScalarMatrix> omegaEqn
322  (
323  fvm::ddt(alpha, rho, omega_)
324  + fvm::div(alphaRhoPhi, omega_)
325  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
326  ==
327  gamma_*alpha()*rho()*G*omega_()/k_()
328  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
329  - fvm::Sp(beta(gradU)*alpha()*rho()*omega_(), omega_)
330  + alpha()*rho()*CDkOmega()
331  + omegaSource()
332  + fvModels.source(alpha, rho, omega_)
333  );
334 
335  omegaEqn.ref().relax();
336  fvConstraints.constrain(omegaEqn.ref());
337  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
338  solve(omegaEqn);
339  fvConstraints.constrain(omega_);
340  boundOmega();
341 
342 
343  // Turbulent kinetic energy equation
345  (
346  fvm::ddt(alpha, rho, k_)
347  + fvm::div(alphaRhoPhi, k_)
348  - fvm::laplacian(alpha*rho*DkEff(), k_)
349  ==
350  alpha()*rho()*G
351  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
352  - fvm::Sp(betaStar_*alpha()*rho()*omega_(), k_)
353  + kSource()
354  + fvModels.source(alpha, rho, k_)
355  );
356 
357  kEqn.ref().relax();
358  fvConstraints.constrain(kEqn.ref());
359  solve(kEqn);
361  bound(k_, this->kMin_);
362  boundOmega();
363 
364  correctNut(gradU);
365 }
366 
367 
368 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
369 
370 } // End namespace RASModels
371 } // End namespace Foam
372 
373 // ************************************************************************* //
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.
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:31
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
void boundOmega()
Bound omega.
Definition: kOmega2006.C:41
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega2006.C:283
virtual tmp< fvScalarMatrix > omegaSource() const
Source term for the omega equation.
Definition: kOmega2006.C:122
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: kOmega2006.C:100
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
Definition: kOmega2006.C:107
kOmega2006(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: kOmega2006.C:139
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega2006.C:262
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
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.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
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, const SuType &Su)
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 > > S(const Pair< tmp< volScalarField::Internal >> &, const VolField< Type > &)
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.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimless
static const Identity< scalar > I
Definition: Identity.H:93
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
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
dimensionedTensor skew(const dimensionedTensor &dt)
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.