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  this->printCoeffs(type);
151  }
152 }
153 
154 
155 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
156 
157 template<class BasicTurbulenceModel>
159 {
161  {
162  Cmu_.readIfPresent(this->coeffDict());
163  beta_.readIfPresent(this->coeffDict());
164  gamma_.readIfPresent(this->coeffDict());
165  alphaK_.readIfPresent(this->coeffDict());
166  alphaOmega_.readIfPresent(this->coeffDict());
167 
168  return true;
169  }
170  else
171  {
172  return false;
173  }
174 }
175 
176 
177 template<class BasicTurbulenceModel>
179 {
180  if (!this->turbulence_)
181  {
182  return;
183  }
184 
185  // Local references
186  const alphaField& alpha = this->alpha_;
187  const rhoField& rho = this->rho_;
188  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
189  const volVectorField& U = this->U_;
190  volScalarField& nut = this->nut_;
191 
193 
195 
196  tmp<volTensorField> tgradU = fvc::grad(U);
198  (
199  this->GName(),
200  nut*(tgradU() && dev(twoSymm(tgradU())))
201  );
202  tgradU.clear();
203 
204  // Update omega and G at the wall
205  omega_.boundaryField().updateCoeffs();
206 
207  // Turbulence specific dissipation rate equation
208  tmp<fvScalarMatrix> omegaEqn
209  (
210  fvm::ddt(alpha, rho, omega_)
211  + fvm::div(alphaRhoPhi, omega_)
212  - fvm::laplacian(alpha*rho*DomegaEff(), omega_)
213  ==
214  gamma_*alpha*rho*G*omega_/k_
215  - fvm::SuSp(((2.0/3.0)*gamma_)*alpha*rho*divU, omega_)
216  - fvm::Sp(beta_*alpha*rho*omega_, omega_)
217  );
218 
219  omegaEqn().relax();
220 
221  omegaEqn().boundaryManipulate(omega_.boundaryField());
222 
223  solve(omegaEqn);
224  bound(omega_, this->omegaMin_);
225 
226 
227  // Turbulent kinetic energy equation
229  (
230  fvm::ddt(alpha, rho, k_)
231  + fvm::div(alphaRhoPhi, k_)
232  - fvm::laplacian(alpha*rho*DkEff(), k_)
233  ==
234  alpha*rho*G
235  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
236  - fvm::Sp(Cmu_*alpha*rho*omega_, k_)
237  );
238 
239  kEqn().relax();
240  solve(kEqn);
241  bound(k_, this->kMin_);
242 
243  correctNut();
244 }
245 
246 
247 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
248 
249 } // End namespace RASModels
250 } // End namespace Foam
251 
252 // ************************************************************************* //
Foam::RASModels::kOmega::read
virtual bool read()
Read RASProperties dictionary.
Definition: kOmega.C:158
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::RASModels::kOmega::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: kOmega.H:106
Foam::IOobject::groupName
static word groupName(Name name, const word &group)
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::RASModels::kOmega::correctNut
virtual void correctNut()
Definition: kOmega.C:40
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
Foam::RASModels::kOmega::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: kOmega.H:107
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::RASModels::kOmega::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kOmega.C:178
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::RASModels::kOmega::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: kOmega.H:108
nut
nut
Definition: createTDFields.H:71
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::RASModels::kOmega::kOmega
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
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
rho
rho
Definition: pEqn.H:3
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
kOmega.H
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104