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-2021 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  this->nut_ = k_/omega_;
44  this->nut_.correctBoundaryConditions();
45  fvConstraints::New(this->mesh_).constrain(this->nut_);
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  return tmp<fvScalarMatrix>
53  (
54  new fvScalarMatrix
55  (
56  k_,
57  dimVolume*this->rho_.dimensions()*k_.dimensions()
58  /dimTime
59  )
60  );
61 }
62 
63 
64 template<class BasicMomentumTransportModel>
66 {
67  return tmp<fvScalarMatrix>
68  (
69  new fvScalarMatrix
70  (
71  omega_,
72  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
73  )
74  );
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
80 template<class BasicMomentumTransportModel>
82 (
83  const alphaField& alpha,
84  const rhoField& rho,
85  const volVectorField& U,
86  const surfaceScalarField& alphaRhoPhi,
87  const surfaceScalarField& phi,
88  const transportModel& transport,
89  const word& type
90 )
91 :
93  (
94  type,
95  alpha,
96  rho,
97  U,
98  alphaRhoPhi,
99  phi,
100  transport
101  ),
102 
103  Cmu_
104  (
106  (
107  "betaStar",
108  this->coeffDict_,
109  0.09
110  )
111  ),
112  beta_
113  (
115  (
116  "beta",
117  this->coeffDict_,
118  0.072
119  )
120  ),
121  gamma_
122  (
124  (
125  "gamma",
126  this->coeffDict_,
127  0.52
128  )
129  ),
130  alphaK_
131  (
133  (
134  "alphaK",
135  this->coeffDict_,
136  0.5
137  )
138  ),
139  alphaOmega_
140  (
142  (
143  "alphaOmega",
144  this->coeffDict_,
145  0.5
146  )
147  ),
148 
149  k_
150  (
151  IOobject
152  (
153  IOobject::groupName("k", alphaRhoPhi.group()),
154  this->runTime_.timeName(),
155  this->mesh_,
158  ),
159  this->mesh_
160  ),
161  omega_
162  (
163  IOobject
164  (
165  IOobject::groupName("omega", alphaRhoPhi.group()),
166  this->runTime_.timeName(),
167  this->mesh_,
170  ),
171  this->mesh_
172  )
173 {
174  bound(k_, this->kMin_);
175  bound(omega_, this->omegaMin_);
176 
177  if (type == typeName)
178  {
179  this->printCoeffs(type);
180  }
181 }
182 
183 
184 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
185 
186 template<class BasicMomentumTransportModel>
188 {
190  {
191  Cmu_.readIfPresent(this->coeffDict());
192  beta_.readIfPresent(this->coeffDict());
193  gamma_.readIfPresent(this->coeffDict());
194  alphaK_.readIfPresent(this->coeffDict());
195  alphaOmega_.readIfPresent(this->coeffDict());
196 
197  return true;
198  }
199  else
200  {
201  return false;
202  }
203 }
204 
205 
206 template<class BasicMomentumTransportModel>
208 {
209  if (!this->turbulence_)
210  {
211  return;
212  }
213 
214  // Local references
215  const alphaField& alpha = this->alpha_;
216  const rhoField& rho = this->rho_;
217  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
218  const volVectorField& U = this->U_;
219  volScalarField& nut = this->nut_;
220  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
222  (
223  Foam::fvConstraints::New(this->mesh_)
224  );
225 
227 
229  (
230  fvc::div(fvc::absolute(this->phi(), U))().v()
231  );
232 
233  tmp<volTensorField> tgradU = fvc::grad(U);
235  (
236  this->GName(),
237  nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v())
238  );
239  tgradU.clear();
240 
241  // Update omega and G at the wall
242  omega_.boundaryFieldRef().updateCoeffs();
243 
244  // Turbulence specific dissipation rate equation
245  tmp<fvScalarMatrix> omegaEqn
246  (
247  fvm::ddt(alpha, rho, omega_)
248  + fvm::div(alphaRhoPhi, omega_)
249  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
250  ==
251  gamma_*alpha()*rho()*G*omega_()/k_()
252  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
253  - fvm::Sp(beta_*alpha()*rho()*omega_(), omega_)
254  + omegaSource()
255  + fvModels.source(alpha, rho, omega_)
256  );
257 
258  omegaEqn.ref().relax();
259  fvConstraints.constrain(omegaEqn.ref());
260  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
261  solve(omegaEqn);
262  fvConstraints.constrain(omega_);
263  bound(omega_, this->omegaMin_);
264 
265 
266  // Turbulent kinetic energy equation
268  (
269  fvm::ddt(alpha, rho, k_)
270  + fvm::div(alphaRhoPhi, k_)
271  - fvm::laplacian(alpha*rho*DkEff(), k_)
272  ==
273  alpha()*rho()*G
274  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
275  - fvm::Sp(Cmu_*alpha()*rho()*omega_(), k_)
276  + kSource()
277  + fvModels.source(alpha, rho, k_)
278  );
279 
280  kEqn.ref().relax();
281  fvConstraints.constrain(kEqn.ref());
282  solve(kEqn);
283  fvConstraints.constrain(k_);
284  bound(k_, this->kMin_);
285 
286  correctNut();
287 }
288 
289 
290 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
291 
292 } // End namespace RASModels
293 } // End namespace Foam
294 
295 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
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
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
Generic dimensioned Type class.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
BasicMomentumTransportModel::rhoField rhoField
Definition: kOmega.H:106
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
const dimensionSet dimTime
bool read(const char *, int32_t &)
Definition: int32IO.C:85
const Internal & v() const
Return a const-reference to the dimensioned internal field.
Definition: volFieldsI.H:31
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
static word groupName(Name name, const word &group)
Finite volume models.
Definition: fvModels.H:60
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
phi
Definition: correctPhi.H:3
BasicMomentumTransportModel::alphaField alphaField
Definition: kOmega.H:105
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:207
kOmega(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: kOmega.C:82
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Bound the given scalar field if it has gone unbounded.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
virtual void correctNut()
Definition: kOmega.C:41
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:187
Foam::fvModels & fvModels
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Finite volume constraints.
Definition: fvConstraints.H:57
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
BasicMomentumTransportModel::transportModel transportModel
Definition: kOmega.H:107
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmega.C:50
const dimensionSet dimVolume
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
scalar nut
Namespace for OpenFOAM.
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmega.C:65