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-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 "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_("betaStar", this->coeffDict(), 0.09),
111  beta_("beta", this->coeffDict(), 0.072),
112  gamma_("gamma", this->coeffDict(), 0.52),
113  alphaK_("alphaK", this->coeffDict(), 0.5),
114  alphaOmega_("alphaOmega", this->coeffDict(), 0.5),
115 
116  k_
117  (
118  IOobject
119  (
120  this->groupName("k"),
121  this->runTime_.name(),
122  this->mesh_,
123  IOobject::MUST_READ,
124  IOobject::AUTO_WRITE
125  ),
126  this->mesh_
127  ),
128  omega_
129  (
130  IOobject
131  (
132  this->groupName("omega"),
133  this->runTime_.name(),
134  this->mesh_,
135  IOobject::MUST_READ,
136  IOobject::AUTO_WRITE
137  ),
138  this->mesh_
139  )
140 {
141  bound(k_, this->kMin_);
142  boundOmega();
143 }
144 
145 
146 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
147 
148 template<class BasicMomentumTransportModel>
150 {
152  {
153  betaStar_.readIfPresent(this->coeffDict());
154  beta_.readIfPresent(this->coeffDict());
155  gamma_.readIfPresent(this->coeffDict());
156  alphaK_.readIfPresent(this->coeffDict());
157  alphaOmega_.readIfPresent(this->coeffDict());
158 
159  return true;
160  }
161  else
162  {
163  return false;
164  }
165 }
166 
167 
168 template<class BasicMomentumTransportModel>
170 {
171  if (!this->turbulence_)
172  {
173  return;
174  }
175 
176  // Local references
177  const alphaField& alpha = this->alpha_;
178  const rhoField& rho = this->rho_;
179  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
180  const volVectorField& U = this->U_;
181  volScalarField& nut = this->nut_;
182  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
184  (
185  Foam::fvConstraints::New(this->mesh_)
186  );
187 
189 
191  (
192  fvc::div(fvc::absolute(this->phi(), U))().v()
193  );
194 
195  tmp<volTensorField> tgradU = fvc::grad(U);
197  (
198  this->GName(),
199  nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v())
200  );
201  tgradU.clear();
202 
203  // Update omega and G at the wall
204  omega_.boundaryFieldRef().updateCoeffs();
205 
206  // Turbulence specific dissipation rate equation
207  tmp<fvScalarMatrix> omegaEqn
208  (
209  fvm::ddt(alpha, rho, omega_)
210  + fvm::div(alphaRhoPhi, omega_)
211  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
212  ==
213  gamma_*alpha()*rho()*G*omega_()/k_()
214  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
215  - fvm::Sp(beta_*alpha()*rho()*omega_(), omega_)
216  + omegaSource()
217  + fvModels.source(alpha, rho, omega_)
218  );
219 
220  omegaEqn.ref().relax();
221  fvConstraints.constrain(omegaEqn.ref());
222  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
223  solve(omegaEqn);
224  fvConstraints.constrain(omega_);
225  boundOmega();
226 
227 
228  // Turbulent kinetic energy equation
230  (
231  fvm::ddt(alpha, rho, k_)
232  + fvm::div(alphaRhoPhi, k_)
233  - fvm::laplacian(alpha*rho*DkEff(), k_)
234  ==
235  alpha()*rho()*G
236  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
237  - fvm::Sp(betaStar_*alpha()*rho()*omega_(), k_)
238  + kSource()
239  + fvModels.source(alpha, rho, k_)
240  );
241 
242  kEqn.ref().relax();
243  fvConstraints.constrain(kEqn.ref());
244  solve(kEqn);
246  bound(k_, this->kMin_);
247  boundOmega();
248 
249  correctNut();
250 }
251 
252 
253 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
254 
255 } // End namespace RASModels
256 } // End namespace Foam
257 
258 // ************************************************************************* //
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:169
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:149
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
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:253
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 > > 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.
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 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
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
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.