kOmega.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | Copyright (C) 2011-2016 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 BasicTurbulenceModel>
41 {
42  this->nut_ = k_/omega_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 
46  BasicTurbulenceModel::correctNut();
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
52 template<class BasicTurbulenceModel>
54 (
55  const alphaField& alpha,
56  const rhoField& rho,
57  const volVectorField& U,
58  const surfaceScalarField& alphaRhoPhi,
59  const surfaceScalarField& phi,
60  const transportModel& transport,
61  const word& propertiesName,
62  const word& type
63 )
64 :
66  (
67  type,
68  alpha,
69  rho,
70  U,
71  alphaRhoPhi,
72  phi,
73  transport,
74  propertiesName
75  ),
76 
77  Cmu_
78  (
80  (
81  "betaStar",
82  this->coeffDict_,
83  0.09
84  )
85  ),
86  beta_
87  (
89  (
90  "beta",
91  this->coeffDict_,
92  0.072
93  )
94  ),
95  gamma_
96  (
98  (
99  "gamma",
100  this->coeffDict_,
101  0.52
102  )
103  ),
104  alphaK_
105  (
107  (
108  "alphaK",
109  this->coeffDict_,
110  0.5
111  )
112  ),
113  alphaOmega_
114  (
116  (
117  "alphaOmega",
118  this->coeffDict_,
119  0.5
120  )
121  ),
122 
123  k_
124  (
125  IOobject
126  (
127  IOobject::groupName("k", U.group()),
128  this->runTime_.timeName(),
129  this->mesh_,
132  ),
133  this->mesh_
134  ),
135  omega_
136  (
137  IOobject
138  (
139  IOobject::groupName("omega", U.group()),
140  this->runTime_.timeName(),
141  this->mesh_,
144  ),
145  this->mesh_
146  )
147 {
148  bound(k_, this->kMin_);
149  bound(omega_, this->omegaMin_);
150 
151  if (type == typeName)
152  {
153  this->printCoeffs(type);
154  }
155 }
156 
157 
158 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
159 
160 template<class BasicTurbulenceModel>
162 {
164  {
165  Cmu_.readIfPresent(this->coeffDict());
166  beta_.readIfPresent(this->coeffDict());
167  gamma_.readIfPresent(this->coeffDict());
168  alphaK_.readIfPresent(this->coeffDict());
169  alphaOmega_.readIfPresent(this->coeffDict());
170 
171  return true;
172  }
173  else
174  {
175  return false;
176  }
177 }
178 
179 
180 template<class BasicTurbulenceModel>
182 {
183  if (!this->turbulence_)
184  {
185  return;
186  }
187 
188  // Local references
189  const alphaField& alpha = this->alpha_;
190  const rhoField& rho = this->rho_;
191  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
192  const volVectorField& U = this->U_;
193  volScalarField& nut = this->nut_;
194  fv::options& fvOptions(fv::options::New(this->mesh_));
195 
197 
199 
200  tmp<volTensorField> tgradU = fvc::grad(U);
202  (
203  this->GName(),
204  nut*(tgradU() && dev(twoSymm(tgradU())))
205  );
206  tgradU.clear();
207 
208  // Update omega and G at the wall
209  omega_.boundaryFieldRef().updateCoeffs();
210 
211  // Turbulence specific dissipation rate equation
212  tmp<fvScalarMatrix> omegaEqn
213  (
214  fvm::ddt(alpha, rho, omega_)
215  + fvm::div(alphaRhoPhi, omega_)
216  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
217  ==
218  gamma_*alpha*rho*G*omega_/k_
219  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha*rho*divU, omega_)
220  - fvm::Sp(beta_*alpha*rho*omega_, omega_)
221  + fvOptions(alpha, rho, omega_)
222  );
223 
224  omegaEqn.ref().relax();
225  fvOptions.constrain(omegaEqn.ref());
226  omegaEqn.ref().boundaryManipulate(omega_.boundaryFieldRef());
227  solve(omegaEqn);
228  fvOptions.correct(omega_);
229  bound(omega_, this->omegaMin_);
230 
231 
232  // Turbulent kinetic energy equation
234  (
235  fvm::ddt(alpha, rho, k_)
236  + fvm::div(alphaRhoPhi, k_)
237  - fvm::laplacian(alpha*rho*DkEff(), k_)
238  ==
239  alpha*rho*G
240  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
241  - fvm::Sp(Cmu_*alpha*rho*omega_, k_)
242  + fvOptions(alpha, rho, k_)
243  );
244 
245  kEqn.ref().relax();
246  fvOptions.constrain(kEqn.ref());
247  solve(kEqn);
248  fvOptions.correct(k_);
249  bound(k_, this->kMin_);
250 
251  correctNut();
252 }
253 
254 
255 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
256 
257 } // End namespace RASModels
258 } // End namespace Foam
259 
260 // ************************************************************************* //
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
U
Definition: pEqn.H:83
BasicTurbulenceModel::transportModel transportModel
Definition: kOmega.H:108
kOmega(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &propertiesName=turbulenceModel::propertiesName, const word &type=typeName)
Construct from components.
Definition: kOmega.C:54
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Generic dimensioned Type class.
BasicTurbulenceModel::alphaField alphaField
Definition: kOmega.H:106
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite-volume options.
Definition: fvOptions.H:52
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fv::options & fvOptions
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
BasicTurbulenceModel::rhoField rhoField
Definition: kOmega.H:107
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:181
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
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(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil)-psi)*pSat, rhoMin);#1"/home/ubuntu/OpenFOAM-4.1/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:72
virtual void correctNut()
Definition: kOmega.C:40
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:161
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:461
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
A class for managing temporary objects.
Definition: PtrList.H:54
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
scalar nut
Namespace for OpenFOAM.