RNGkEpsilon.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 "RNGkEpsilon.H"
27 #include "bound.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace RASModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicTurbulenceModel>
40 {
41  this->nut_ = Cmu_*sqr(k_)/epsilon_;
42  this->nut_.correctBoundaryConditions();
43 
44  BasicTurbulenceModel::correctNut();
45 }
46 
47 
48 template<class BasicTurbulenceModel>
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 BasicTurbulenceModel>
65 {
66  return tmp<fvScalarMatrix>
67  (
68  new fvScalarMatrix
69  (
70  epsilon_,
71  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
72  /dimTime
73  )
74  );
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
80 template<class BasicTurbulenceModel>
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& propertiesName,
90  const word& type
91 )
92 :
94  (
95  type,
96  alpha,
97  rho,
98  U,
99  alphaRhoPhi,
100  phi,
101  transport,
102  propertiesName
103  ),
104 
105  Cmu_
106  (
108  (
109  "Cmu",
110  this->coeffDict_,
111  0.0845
112  )
113  ),
114  C1_
115  (
117  (
118  "C1",
119  this->coeffDict_,
120  1.42
121  )
122  ),
123  C2_
124  (
126  (
127  "C2",
128  this->coeffDict_,
129  1.68
130  )
131  ),
132  C3_
133  (
135  (
136  "C3",
137  this->coeffDict_,
138  -0.33
139  )
140  ),
141  sigmak_
142  (
144  (
145  "sigmak",
146  this->coeffDict_,
147  0.71942
148  )
149  ),
150  sigmaEps_
151  (
153  (
154  "sigmaEps",
155  this->coeffDict_,
156  0.71942
157  )
158  ),
159  eta0_
160  (
162  (
163  "eta0",
164  this->coeffDict_,
165  4.38
166  )
167  ),
168  beta_
169  (
171  (
172  "beta",
173  this->coeffDict_,
174  0.012
175  )
176  ),
177 
178  k_
179  (
180  IOobject
181  (
182  IOobject::groupName("k", U.group()),
183  this->runTime_.timeName(),
184  this->mesh_,
187  ),
188  this->mesh_
189  ),
190  epsilon_
191  (
192  IOobject
193  (
194  IOobject::groupName("epsilon", U.group()),
195  this->runTime_.timeName(),
196  this->mesh_,
199  ),
200  this->mesh_
201  )
202 {
203  bound(k_, this->kMin_);
204  bound(epsilon_, this->epsilonMin_);
205 
206  if (type == typeName)
207  {
208  this->printCoeffs(type);
209  }
210 }
211 
212 
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
214 
215 template<class BasicTurbulenceModel>
217 {
219  {
220  Cmu_.readIfPresent(this->coeffDict());
221  C1_.readIfPresent(this->coeffDict());
222  C2_.readIfPresent(this->coeffDict());
223  C3_.readIfPresent(this->coeffDict());
224  sigmak_.readIfPresent(this->coeffDict());
225  sigmaEps_.readIfPresent(this->coeffDict());
226  eta0_.readIfPresent(this->coeffDict());
227  beta_.readIfPresent(this->coeffDict());
228 
229  return true;
230  }
231  else
232  {
233  return false;
234  }
235 }
236 
237 
238 template<class BasicTurbulenceModel>
240 {
241  if (!this->turbulence_)
242  {
243  return;
244  }
245 
246  // Local references
247  const alphaField& alpha = this->alpha_;
248  const rhoField& rho = this->rho_;
249  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
250  const volVectorField& U = this->U_;
251  volScalarField& nut = this->nut_;
252 
254 
256 
257  tmp<volTensorField> tgradU = fvc::grad(U);
258  volScalarField S2((tgradU() && dev(twoSymm(tgradU()))));
259  tgradU.clear();
260 
261  volScalarField G(this->GName(), nut*S2);
262 
263  volScalarField eta(sqrt(mag(S2))*k_/epsilon_);
264  volScalarField eta3(eta*sqr(eta));
265 
267  (
268  ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
269  );
270 
271  // Update epsilon and G at the wall
272  epsilon_.boundaryField().updateCoeffs();
273 
274  // Dissipation equation
275  tmp<fvScalarMatrix> epsEqn
276  (
277  fvm::ddt(alpha, rho, epsilon_)
278  + fvm::div(alphaRhoPhi, epsilon_)
279  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
280  ==
281  (C1_ - R)*alpha*rho*G*epsilon_/k_
282  - fvm::SuSp(((2.0/3.0)*C1_ + C3_)*alpha*rho*divU, epsilon_)
283  - fvm::Sp(C2_*alpha*rho*epsilon_/k_, epsilon_)
284  + epsilonSource()
285  );
286 
287  epsEqn().relax();
288 
289  epsEqn().boundaryManipulate(epsilon_.boundaryField());
290 
291  solve(epsEqn);
292  bound(epsilon_, this->epsilonMin_);
293 
294 
295  // Turbulent kinetic energy equation
296 
298  (
299  fvm::ddt(alpha, rho, k_)
300  + fvm::div(alphaRhoPhi, k_)
301  - fvm::laplacian(alpha*rho*DkEff(), k_)
302  ==
303  alpha*rho*G
304  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
305  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
306  + kSource()
307  );
308 
309  kEqn().relax();
310  solve(kEqn);
311  bound(k_, this->kMin_);
312 
313  correctNut();
314 }
315 
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 } // End namespace RASModels
320 } // End namespace Foam
321 
322 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
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::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
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::RASModels::RNGkEpsilon::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:49
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::RNGkEpsilon::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:64
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
R
#define R(A, B, C, D, E, F, K, M)
Foam::RASModels::RNGkEpsilon::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:239
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::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
Foam::RASModels::RNGkEpsilon::correctNut
virtual void correctNut()
Definition: RNGkEpsilon.C:39
rho
rho
Definition: pEqn.H:3
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
Foam::RASModels::RNGkEpsilon::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:216
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
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
RNGkEpsilon.H
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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::RASModels::RNGkEpsilon::RNGkEpsilon
RNGkEpsilon(const RNGkEpsilon &)
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104