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  correctNut();
209  this->printCoeffs(type);
210  }
211 }
212 
213 
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 
216 template<class BasicTurbulenceModel>
218 {
220  {
221  Cmu_.readIfPresent(this->coeffDict());
222  C1_.readIfPresent(this->coeffDict());
223  C2_.readIfPresent(this->coeffDict());
224  C3_.readIfPresent(this->coeffDict());
225  sigmak_.readIfPresent(this->coeffDict());
226  sigmaEps_.readIfPresent(this->coeffDict());
227  eta0_.readIfPresent(this->coeffDict());
228  beta_.readIfPresent(this->coeffDict());
229 
230  return true;
231  }
232  else
233  {
234  return false;
235  }
236 }
237 
238 
239 template<class BasicTurbulenceModel>
241 {
242  if (!this->turbulence_)
243  {
244  return;
245  }
246 
247  // Local references
248  const alphaField& alpha = this->alpha_;
249  const rhoField& rho = this->rho_;
250  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
251  const volVectorField& U = this->U_;
252  volScalarField& nut = this->nut_;
253 
255 
257 
258  tmp<volTensorField> tgradU = fvc::grad(U);
259  volScalarField S2((tgradU() && dev(twoSymm(tgradU()))));
260  tgradU.clear();
261 
262  volScalarField G(this->GName(), nut*S2);
263 
264  volScalarField eta(sqrt(mag(S2))*k_/epsilon_);
265  volScalarField eta3(eta*sqr(eta));
266 
268  (
269  ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
270  );
271 
272  // Update epsilon and G at the wall
273  epsilon_.boundaryField().updateCoeffs();
274 
275  // Dissipation equation
276  tmp<fvScalarMatrix> epsEqn
277  (
278  fvm::ddt(alpha, rho, epsilon_)
279  + fvm::div(alphaRhoPhi, epsilon_)
280  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
281  ==
282  (C1_ - R)*alpha*rho*G*epsilon_/k_
283  - fvm::SuSp(((2.0/3.0)*C1_ + C3_)*alpha*rho*divU, epsilon_)
284  - fvm::Sp(C2_*alpha*rho*epsilon_/k_, epsilon_)
285  + epsilonSource()
286  );
287 
288  epsEqn().relax();
289 
290  epsEqn().boundaryManipulate(epsilon_.boundaryField());
291 
292  solve(epsEqn);
293  bound(epsilon_, this->epsilonMin_);
294 
295 
296  // Turbulent kinetic energy equation
297 
299  (
300  fvm::ddt(alpha, rho, k_)
301  + fvm::div(alphaRhoPhi, k_)
302  - fvm::laplacian(alpha*rho*DkEff(), k_)
303  ==
304  alpha*rho*G
305  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
306  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
307  + kSource()
308  );
309 
310  kEqn().relax();
311  solve(kEqn);
312  bound(k_, this->kMin_);
313 
314  correctNut();
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 } // End namespace RASModels
321 } // End namespace Foam
322 
323 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
Renormalization group k-epsilon turbulence model for incompressible and compressible flows...
Definition: RNGkEpsilon.H:85
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:217
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
dimensioned< scalar > mag(const dimensioned< Type > &)
Bound the given scalar field if it has gone unbounded.
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:174
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
#define R(A, B, C, D, E, F, K, M)
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
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:64
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)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
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 > &)
scalar nut
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:49
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
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)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:240
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
bool read(const char *, int32_t &)
Definition: int32IO.C:87
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99
A class for managing temporary objects.
Definition: PtrList.H:118
const dimensionedScalar G
Newtonian constant of gravitation.
U
Definition: pEqn.H:82