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-2024 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_("betaStar", this->coeffDict(), 0.09),
161  beta0_("beta0", this->coeffDict(), 0.0708),
162  gamma_("gamma", this->coeffDict(), 0.52),
163  Clim_("Clim", this->coeffDict(), 0.875),
164  sigmaDo_("sigmaDo", this->coeffDict(), 0.125),
165  alphaK_("alphaK", this->coeffDict(), 0.6),
166  alphaOmega_("alphaOmega", this->coeffDict(), 0.5),
167 
168  k_
169  (
170  IOobject
171  (
172  this->groupName("k"),
173  this->runTime_.name(),
174  this->mesh_,
175  IOobject::MUST_READ,
176  IOobject::AUTO_WRITE
177  ),
178  this->mesh_
179  ),
180  omega_
181  (
182  IOobject
183  (
184  this->groupName("omega"),
185  this->runTime_.name(),
186  this->mesh_,
187  IOobject::MUST_READ,
188  IOobject::AUTO_WRITE
189  ),
190  this->mesh_
191  )
192 {
193  bound(k_, this->kMin_);
194  boundOmega();
195 }
196 
197 
198 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
199 
200 template<class BasicMomentumTransportModel>
202 {
204  {
205  betaStar_.readIfPresent(this->coeffDict());
206  beta0_.readIfPresent(this->coeffDict());
207  gamma_.readIfPresent(this->coeffDict());
208  sigmaDo_.readIfPresent(this->coeffDict());
209  alphaK_.readIfPresent(this->coeffDict());
210  alphaOmega_.readIfPresent(this->coeffDict());
211 
212  return true;
213  }
214  else
215  {
216  return false;
217  }
218 }
219 
220 
221 template<class BasicMomentumTransportModel>
223 {
224  if (!this->turbulence_)
225  {
226  return;
227  }
228 
229  // Local references
230  const alphaField& alpha = this->alpha_;
231  const rhoField& rho = this->rho_;
232  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
233  const volVectorField& U = this->U_;
234  volScalarField& nut = this->nut_;
235  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
237  (
238  Foam::fvConstraints::New(this->mesh_)
239  );
240 
242 
244  (
245  typedName("divU"),
246  fvc::div(fvc::absolute(this->phi(), U))().v()
247  );
248 
249  const volTensorField gradU(fvc::grad(U));
251  (
252  this->GName(),
253  nut.v()*(dev(twoSymm(gradU.v())) && gradU.v())
254  );
255 
256  // Update omega and G at the wall
257  omega_.boundaryFieldRef().updateCoeffs();
258 
259  // Turbulence specific dissipation rate equation
260  tmp<fvScalarMatrix> omegaEqn
261  (
262  fvm::ddt(alpha, rho, omega_)
263  + fvm::div(alphaRhoPhi, omega_)
264  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
265  ==
266  gamma_*alpha()*rho()*G*omega_()/k_()
267  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
268  - fvm::Sp(beta(gradU)*alpha()*rho()*omega_(), omega_)
269  + alpha()*rho()*CDkOmega()
270  + omegaSource()
271  + fvModels.source(alpha, rho, omega_)
272  );
273 
274  omegaEqn.ref().relax();
275  fvConstraints.constrain(omegaEqn.ref());
276  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
277  solve(omegaEqn);
278  fvConstraints.constrain(omega_);
279  boundOmega();
280 
281 
282  // Turbulent kinetic energy equation
284  (
285  fvm::ddt(alpha, rho, k_)
286  + fvm::div(alphaRhoPhi, k_)
287  - fvm::laplacian(alpha*rho*DkEff(), k_)
288  ==
289  alpha()*rho()*G
290  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
291  - fvm::Sp(betaStar_*alpha()*rho()*omega_(), k_)
292  + kSource()
293  + fvModels.source(alpha, rho, k_)
294  );
295 
296  kEqn.ref().relax();
297  fvConstraints.constrain(kEqn.ref());
298  solve(kEqn);
300  bound(k_, this->kMin_);
301  boundOmega();
302 
303  correctNut(gradU);
304 }
305 
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
309 } // End namespace RASModels
310 } // End namespace Foam
311 
312 // ************************************************************************* //
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:29
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:222
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:201
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
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:197
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 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.
void skew(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimless
static const Identity< scalar > I
Definition: Identity.H:93
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
void symm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
const dimensionSet dimVolume
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
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
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void tr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< tensor > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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.