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-2020 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 "fvOptions.H"
28 #include "bound.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace RASModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicMomentumTransportModel>
41 {
42  this->nut_ = k_/omega_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 }
46 
47 
48 template<class BasicMomentumTransportModel>
50 {
51  return tmp<fvScalarMatrix>
52  (
53  new fvScalarMatrix
54  (
55  k_,
56  dimVolume*this->rho_.dimensions()*k_.dimensions()
57  /dimTime
58  )
59  );
60 }
61 
62 
63 template<class BasicMomentumTransportModel>
65 {
66  return tmp<fvScalarMatrix>
67  (
68  new fvScalarMatrix
69  (
70  omega_,
71  dimVolume*this->rho_.dimensions()*omega_.dimensions()/dimTime
72  )
73  );
74 }
75 
76 
77 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
78 
79 template<class BasicMomentumTransportModel>
81 (
82  const alphaField& alpha,
83  const rhoField& rho,
84  const volVectorField& U,
85  const surfaceScalarField& alphaRhoPhi,
86  const surfaceScalarField& phi,
87  const transportModel& transport,
88  const word& type
89 )
90 :
92  (
93  type,
94  alpha,
95  rho,
96  U,
97  alphaRhoPhi,
98  phi,
99  transport
100  ),
101 
102  Cmu_
103  (
105  (
106  "betaStar",
107  this->coeffDict_,
108  0.09
109  )
110  ),
111  beta_
112  (
114  (
115  "beta",
116  this->coeffDict_,
117  0.072
118  )
119  ),
120  gamma_
121  (
123  (
124  "gamma",
125  this->coeffDict_,
126  0.52
127  )
128  ),
129  alphaK_
130  (
132  (
133  "alphaK",
134  this->coeffDict_,
135  0.5
136  )
137  ),
138  alphaOmega_
139  (
141  (
142  "alphaOmega",
143  this->coeffDict_,
144  0.5
145  )
146  ),
147 
148  k_
149  (
150  IOobject
151  (
152  IOobject::groupName("k", alphaRhoPhi.group()),
153  this->runTime_.timeName(),
154  this->mesh_,
157  ),
158  this->mesh_
159  ),
160  omega_
161  (
162  IOobject
163  (
164  IOobject::groupName("omega", alphaRhoPhi.group()),
165  this->runTime_.timeName(),
166  this->mesh_,
169  ),
170  this->mesh_
171  )
172 {
173  bound(k_, this->kMin_);
174  bound(omega_, this->omegaMin_);
175 
176  if (type == typeName)
177  {
178  this->printCoeffs(type);
179  }
180 }
181 
182 
183 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
184 
185 template<class BasicMomentumTransportModel>
187 {
189  {
190  Cmu_.readIfPresent(this->coeffDict());
191  beta_.readIfPresent(this->coeffDict());
192  gamma_.readIfPresent(this->coeffDict());
193  alphaK_.readIfPresent(this->coeffDict());
194  alphaOmega_.readIfPresent(this->coeffDict());
195 
196  return true;
197  }
198  else
199  {
200  return false;
201  }
202 }
203 
204 
205 template<class BasicMomentumTransportModel>
207 {
208  if (!this->turbulence_)
209  {
210  return;
211  }
212 
213  // Local references
214  const alphaField& alpha = this->alpha_;
215  const rhoField& rho = this->rho_;
216  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
217  const volVectorField& U = this->U_;
218  volScalarField& nut = this->nut_;
219  fv::options& fvOptions(fv::options::New(this->mesh_));
220 
222 
224  (
225  fvc::div(fvc::absolute(this->phi(), U))().v()
226  );
227 
228  tmp<volTensorField> tgradU = fvc::grad(U);
230  (
231  this->GName(),
232  nut.v()*(dev(twoSymm(tgradU().v())) && tgradU().v())
233  );
234  tgradU.clear();
235 
236  // Update omega and G at the wall
237  omega_.boundaryFieldRef().updateCoeffs();
238 
239  // Turbulence specific dissipation rate equation
240  tmp<fvScalarMatrix> omegaEqn
241  (
242  fvm::ddt(alpha, rho, omega_)
243  + fvm::div(alphaRhoPhi, omega_)
244  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
245  ==
246  gamma_*alpha()*rho()*G*omega_()/k_()
247  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha()*rho()*divU, omega_)
248  - fvm::Sp(beta_*alpha()*rho()*omega_(), omega_)
249  + omegaSource()
250  + fvOptions(alpha, rho, omega_)
251  );
252 
253  omegaEqn.ref().relax();
254  fvOptions.constrain(omegaEqn.ref());
255  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
256  solve(omegaEqn);
257  fvOptions.correct(omega_);
258  bound(omega_, this->omegaMin_);
259 
260 
261  // Turbulent kinetic energy equation
263  (
264  fvm::ddt(alpha, rho, k_)
265  + fvm::div(alphaRhoPhi, k_)
266  - fvm::laplacian(alpha*rho*DkEff(), k_)
267  ==
268  alpha()*rho()*G
269  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
270  - fvm::Sp(Cmu_*alpha()*rho()*omega_(), k_)
271  + kSource()
272  + fvOptions(alpha, rho, k_)
273  );
274 
275  kEqn.ref().relax();
276  fvOptions.constrain(kEqn.ref());
277  solve(kEqn);
278  fvOptions.correct(k_);
279  bound(k_, this->kMin_);
280 
281  correctNut();
282 }
283 
284 
285 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
286 
287 } // End namespace RASModels
288 } // End namespace Foam
289 
290 // ************************************************************************* //
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:230
fv::options & fvOptions
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Generic dimensioned Type class.
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
Finite-volume options.
Definition: fvOptions.H:52
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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)
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
BasicMomentumTransportModel::alphaField alphaField
Definition: kOmega.H:105
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:206
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
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:81
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.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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
virtual void correctNut()
Definition: kOmega.C:40
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:186
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
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...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual tmp< fvScalarMatrix > kSource() const
Definition: kOmega.C:49
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
scalar nut
Namespace for OpenFOAM.
virtual tmp< fvScalarMatrix > omegaSource() const
Definition: kOmega.C:64