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-2015 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 "bound.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace RASModels
34 {
35 
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicTurbulenceModel>
41 {
42  this->nut_ = k_/omega_;
43  this->nut_.correctBoundaryConditions();
44 }
45 
46 
47 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
48 
49 template<class BasicTurbulenceModel>
51 (
52  const alphaField& alpha,
53  const rhoField& rho,
54  const volVectorField& U,
55  const surfaceScalarField& alphaRhoPhi,
56  const surfaceScalarField& phi,
57  const transportModel& transport,
58  const word& propertiesName,
59  const word& type
60 )
61 :
63  (
64  type,
65  alpha,
66  rho,
67  U,
68  alphaRhoPhi,
69  phi,
70  transport,
71  propertiesName
72  ),
73 
74  Cmu_
75  (
77  (
78  "betaStar",
79  this->coeffDict_,
80  0.09
81  )
82  ),
83  beta_
84  (
86  (
87  "beta",
88  this->coeffDict_,
89  0.072
90  )
91  ),
92  gamma_
93  (
95  (
96  "gamma",
97  this->coeffDict_,
98  0.52
99  )
100  ),
101  alphaK_
102  (
104  (
105  "alphaK",
106  this->coeffDict_,
107  0.5
108  )
109  ),
110  alphaOmega_
111  (
113  (
114  "alphaOmega",
115  this->coeffDict_,
116  0.5
117  )
118  ),
119 
120  k_
121  (
122  IOobject
123  (
124  IOobject::groupName("k", U.group()),
125  this->runTime_.timeName(),
126  this->mesh_,
129  ),
130  this->mesh_
131  ),
132  omega_
133  (
134  IOobject
135  (
136  IOobject::groupName("omega", U.group()),
137  this->runTime_.timeName(),
138  this->mesh_,
141  ),
142  this->mesh_
143  )
144 {
145  bound(k_, this->kMin_);
146  bound(omega_, this->omegaMin_);
147 
148  if (type == typeName)
149  {
150  correctNut();
151  this->printCoeffs(type);
152  }
153 }
154 
155 
156 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
157 
158 template<class BasicTurbulenceModel>
160 {
162  {
163  Cmu_.readIfPresent(this->coeffDict());
164  beta_.readIfPresent(this->coeffDict());
165  gamma_.readIfPresent(this->coeffDict());
166  alphaK_.readIfPresent(this->coeffDict());
167  alphaOmega_.readIfPresent(this->coeffDict());
168 
169  return true;
170  }
171  else
172  {
173  return false;
174  }
175 }
176 
177 
178 template<class BasicTurbulenceModel>
180 {
181  if (!this->turbulence_)
182  {
183  return;
184  }
185 
186  // Local references
187  const alphaField& alpha = this->alpha_;
188  const rhoField& rho = this->rho_;
189  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
190  const volVectorField& U = this->U_;
191  volScalarField& nut = this->nut_;
192 
194 
196 
197  tmp<volTensorField> tgradU = fvc::grad(U);
199  (
200  this->GName(),
201  nut*(tgradU() && dev(twoSymm(tgradU())))
202  );
203  tgradU.clear();
204 
205  // Update omega and G at the wall
206  omega_.boundaryField().updateCoeffs();
207 
208  // Turbulence specific dissipation rate equation
209  tmp<fvScalarMatrix> omegaEqn
210  (
211  fvm::ddt(alpha, rho, omega_)
212  + fvm::div(alphaRhoPhi, omega_)
213  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
214  ==
215  gamma_*alpha*rho*G*omega_/k_
216  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha*rho*divU, omega_)
217  - fvm::Sp(beta_*alpha*rho*omega_, omega_)
218  );
219 
220  omegaEqn().relax();
221 
222  omegaEqn().boundaryManipulate(omega_.boundaryField());
223 
224  solve(omegaEqn);
225  bound(omega_, this->omegaMin_);
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(Cmu_*alpha*rho*omega_, k_)
238  );
239 
240  kEqn().relax();
241  solve(kEqn);
242  bound(k_, this->kMin_);
243 
244  correctNut();
245 }
246 
247 
248 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
249 
250 } // End namespace RASModels
251 } // End namespace Foam
252 
253 // ************************************************************************* //
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
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.
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:174
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:159
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/dm2/henry/OpenFOAM/OpenFOAM-3.0.x/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:74
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
static word groupName(Name name, const word &group)
BasicTurbulenceModel::rhoField rhoField
Definition: kOmega.H:107
Generic dimensioned Type class.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
BasicTurbulenceModel::transportModel transportModel
Definition: kOmega.H:108
scalar nut
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
BasicTurbulenceModel::alphaField alphaField
Definition: kOmega.H:106
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
surfaceScalarField & phi
solverPerformance solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:87
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:179
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
virtual void correctNut()
Definition: kOmega.C:40
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:51
A class for managing temporary objects.
Definition: PtrList.H:118
const dimensionedScalar G
Newtonian constant of gravitation.
U
Definition: pEqn.H:82