kOmega.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-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 "kOmega.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  this->nut_ = k_/omega_;
51  this->nut_.correctBoundaryConditions();
52  fvConstraints::New(this->mesh_).constrain(this->nut_);
53 }
54 
55 
56 template<class BasicMomentumTransportModel>
58 {
59  return tmp<fvScalarMatrix>
60  (
61  new fvScalarMatrix
62  (
63  k_,
64  dimVolume*this->rho_.dimensions()*k_.dimensions()
65  /dimTime
66  )
67  );
68 }
69 
70 
71 template<class BasicMomentumTransportModel>
73 {
74  return tmp<fvScalarMatrix>
75  (
76  new fvScalarMatrix
77  (
78  omega_,
79  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
80  )
81  );
82 }
83 
84 
85 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
86 
87 template<class BasicMomentumTransportModel>
89 (
90  const alphaField& alpha,
91  const rhoField& rho,
92  const volVectorField& U,
93  const surfaceScalarField& alphaRhoPhi,
94  const surfaceScalarField& phi,
95  const viscosity& viscosity,
96  const word& type
97 )
98 :
99  eddyViscosity<RASModel<BasicMomentumTransportModel>>
100  (
101  type,
102  alpha,
103  rho,
104  U,
105  alphaRhoPhi,
106  phi,
107  viscosity
108  ),
109 
110  betaStar_
111  (
112  dimensioned<scalar>::lookupOrAddToDict
113  (
114  "betaStar",
115  this->coeffDict_,
116  0.09
117  )
118  ),
119  beta_
120  (
121  dimensioned<scalar>::lookupOrAddToDict
122  (
123  "beta",
124  this->coeffDict_,
125  0.072
126  )
127  ),
128  gamma_
129  (
130  dimensioned<scalar>::lookupOrAddToDict
131  (
132  "gamma",
133  this->coeffDict_,
134  0.52
135  )
136  ),
137  alphaK_
138  (
139  dimensioned<scalar>::lookupOrAddToDict
140  (
141  "alphaK",
142  this->coeffDict_,
143  0.5
144  )
145  ),
146  alphaOmega_
147  (
148  dimensioned<scalar>::lookupOrAddToDict
149  (
150  "alphaOmega",
151  this->coeffDict_,
152  0.5
153  )
154  ),
155 
156  k_
157  (
158  IOobject
159  (
160  this->groupName("k"),
161  this->runTime_.name(),
162  this->mesh_,
163  IOobject::MUST_READ,
164  IOobject::AUTO_WRITE
165  ),
166  this->mesh_
167  ),
168  omega_
169  (
170  IOobject
171  (
172  this->groupName("omega"),
173  this->runTime_.name(),
174  this->mesh_,
175  IOobject::MUST_READ,
176  IOobject::AUTO_WRITE
177  ),
178  this->mesh_
179  )
180 {
181  bound(k_, this->kMin_);
182  boundOmega();
183 
184  if (type == typeName)
185  {
186  this->printCoeffs(type);
187  }
188 }
189 
190 
191 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
192 
193 template<class BasicMomentumTransportModel>
195 {
197  {
198  betaStar_.readIfPresent(this->coeffDict());
199  beta_.readIfPresent(this->coeffDict());
200  gamma_.readIfPresent(this->coeffDict());
201  alphaK_.readIfPresent(this->coeffDict());
202  alphaOmega_.readIfPresent(this->coeffDict());
203 
204  return true;
205  }
206  else
207  {
208  return false;
209  }
210 }
211 
212 
213 template<class BasicMomentumTransportModel>
215 {
216  if (!this->turbulence_)
217  {
218  return;
219  }
220 
221  // Local references
222  const alphaField& alpha = this->alpha_;
223  const rhoField& rho = this->rho_;
224  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
225  const volVectorField& U = this->U_;
226  volScalarField& nut = this->nut_;
227  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
229  (
230  Foam::fvConstraints::New(this->mesh_)
231  );
232 
234 
236  (
237  fvc::div(fvc::absolute(this->phi(), U))().v()
238  );
239 
240  tmp<volTensorField> tgradU = fvc::grad(U);
242  (
243  this->GName(),
244  nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v())
245  );
246  tgradU.clear();
247 
248  // Update omega and G at the wall
249  omega_.boundaryFieldRef().updateCoeffs();
250 
251  // Turbulence specific dissipation rate equation
252  tmp<fvScalarMatrix> omegaEqn
253  (
254  fvm::ddt(alpha, rho, omega_)
255  + fvm::div(alphaRhoPhi, omega_)
256  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
257  ==
258  gamma_*alpha()*rho()*G*omega_()/k_()
259  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
260  - fvm::Sp(beta_*alpha()*rho()*omega_(), omega_)
261  + omegaSource()
262  + fvModels.source(alpha, rho, omega_)
263  );
264 
265  omegaEqn.ref().relax();
266  fvConstraints.constrain(omegaEqn.ref());
267  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
268  solve(omegaEqn);
269  fvConstraints.constrain(omega_);
270  boundOmega();
271 
272 
273  // Turbulent kinetic energy equation
275  (
276  fvm::ddt(alpha, rho, k_)
277  + fvm::div(alphaRhoPhi, k_)
278  - fvm::laplacian(alpha*rho*DkEff(), k_)
279  ==
280  alpha()*rho()*G
281  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
282  - fvm::Sp(betaStar_*alpha()*rho()*omega_(), k_)
283  + kSource()
284  + fvModels.source(alpha, rho, k_)
285  );
286 
287  kEqn.ref().relax();
288  fvConstraints.constrain(kEqn.ref());
289  solve(kEqn);
291  bound(k_, this->kMin_);
292  boundOmega();
293 
294  correctNut();
295 }
296 
297 
298 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
299 
300 } // End namespace RASModels
301 } // End namespace Foam
302 
303 // ************************************************************************* //
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.
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: kOmega.C:41
volScalarField k_
Definition: kOmega.H:92
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:214
kOmega(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: kOmega.C:89
virtual tmp< fvScalarMatrix > omegaSource() const
Source term for the omega equation.
Definition: kOmega.C:72
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: kOmega.C:48
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
Definition: kOmega.C:57
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:194
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
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
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 > > 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)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimTime
const dimensionSet dimVolume
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
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.