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-2017 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 "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_ = Cmu_*sqr(k_)/epsilon_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 
46  BasicTurbulenceModel::correctNut();
47 }
48 
49 
50 template<class BasicTurbulenceModel>
52 {
53  return tmp<fvScalarMatrix>
54  (
55  new fvScalarMatrix
56  (
57  k_,
58  dimVolume*this->rho_.dimensions()*k_.dimensions()
59  /dimTime
60  )
61  );
62 }
63 
64 
65 template<class BasicTurbulenceModel>
67 {
68  return tmp<fvScalarMatrix>
69  (
70  new fvScalarMatrix
71  (
72  epsilon_,
73  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
74  /dimTime
75  )
76  );
77 }
78 
79 
80 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
81 
82 template<class BasicTurbulenceModel>
84 (
85  const alphaField& alpha,
86  const rhoField& rho,
87  const volVectorField& U,
88  const surfaceScalarField& alphaRhoPhi,
89  const surfaceScalarField& phi,
90  const transportModel& transport,
91  const word& propertiesName,
92  const word& type
93 )
94 :
96  (
97  type,
98  alpha,
99  rho,
100  U,
101  alphaRhoPhi,
102  phi,
103  transport,
104  propertiesName
105  ),
106 
107  Cmu_
108  (
110  (
111  "Cmu",
112  this->coeffDict_,
113  0.0845
114  )
115  ),
116  C1_
117  (
119  (
120  "C1",
121  this->coeffDict_,
122  1.42
123  )
124  ),
125  C2_
126  (
128  (
129  "C2",
130  this->coeffDict_,
131  1.68
132  )
133  ),
134  C3_
135  (
137  (
138  "C3",
139  this->coeffDict_,
140  0
141  )
142  ),
143  sigmak_
144  (
146  (
147  "sigmak",
148  this->coeffDict_,
149  0.71942
150  )
151  ),
152  sigmaEps_
153  (
155  (
156  "sigmaEps",
157  this->coeffDict_,
158  0.71942
159  )
160  ),
161  eta0_
162  (
164  (
165  "eta0",
166  this->coeffDict_,
167  4.38
168  )
169  ),
170  beta_
171  (
173  (
174  "beta",
175  this->coeffDict_,
176  0.012
177  )
178  ),
179 
180  k_
181  (
182  IOobject
183  (
184  IOobject::groupName("k", U.group()),
185  this->runTime_.timeName(),
186  this->mesh_,
189  ),
190  this->mesh_
191  ),
192  epsilon_
193  (
194  IOobject
195  (
196  IOobject::groupName("epsilon", U.group()),
197  this->runTime_.timeName(),
198  this->mesh_,
201  ),
202  this->mesh_
203  )
204 {
205  bound(k_, this->kMin_);
206  bound(epsilon_, this->epsilonMin_);
207 
208  if (type == typeName)
209  {
210  this->printCoeffs(type);
211  }
212 }
213 
214 
215 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
216 
217 template<class BasicTurbulenceModel>
219 {
221  {
222  Cmu_.readIfPresent(this->coeffDict());
223  C1_.readIfPresent(this->coeffDict());
224  C2_.readIfPresent(this->coeffDict());
225  C3_.readIfPresent(this->coeffDict());
226  sigmak_.readIfPresent(this->coeffDict());
227  sigmaEps_.readIfPresent(this->coeffDict());
228  eta0_.readIfPresent(this->coeffDict());
229  beta_.readIfPresent(this->coeffDict());
230 
231  return true;
232  }
233  else
234  {
235  return false;
236  }
237 }
238 
239 
240 template<class BasicTurbulenceModel>
242 {
243  if (!this->turbulence_)
244  {
245  return;
246  }
247 
248  // Local references
249  const alphaField& alpha = this->alpha_;
250  const rhoField& rho = this->rho_;
251  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
252  const volVectorField& U = this->U_;
253  volScalarField& nut = this->nut_;
254  fv::options& fvOptions(fv::options::New(this->mesh_));
255 
257 
259 
260  tmp<volTensorField> tgradU = fvc::grad(U);
261  volScalarField S2((tgradU() && dev(twoSymm(tgradU()))));
262  tgradU.clear();
263 
264  volScalarField G(this->GName(), nut*S2);
265 
266  volScalarField eta(sqrt(mag(S2))*k_/epsilon_);
267  volScalarField eta3(eta*sqr(eta));
268 
270  (
271  ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
272  );
273 
274  // Update epsilon and G at the wall
275  epsilon_.boundaryFieldRef().updateCoeffs();
276 
277  // Dissipation equation
278  tmp<fvScalarMatrix> epsEqn
279  (
280  fvm::ddt(alpha, rho, epsilon_)
281  + fvm::div(alphaRhoPhi, epsilon_)
282  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
283  ==
284  (C1_ - R)*alpha*rho*G*epsilon_/k_
285  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha*rho*divU, epsilon_)
286  - fvm::Sp(C2_*alpha*rho*epsilon_/k_, epsilon_)
287  + epsilonSource()
288  + fvOptions(alpha, rho, epsilon_)
289  );
290 
291  epsEqn.ref().relax();
292  fvOptions.constrain(epsEqn.ref());
293  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
294  solve(epsEqn);
295  fvOptions.correct(epsilon_);
296  bound(epsilon_, this->epsilonMin_);
297 
298 
299  // Turbulent kinetic energy equation
300 
302  (
303  fvm::ddt(alpha, rho, k_)
304  + fvm::div(alphaRhoPhi, k_)
305  - fvm::laplacian(alpha*rho*DkEff(), k_)
306  ==
307  alpha*rho*G
308  - fvm::SuSp((2.0/3.0)*alpha*rho*divU, k_)
309  - fvm::Sp(alpha*rho*epsilon_/k_, k_)
310  + kSource()
311  + fvOptions(alpha, rho, k_)
312  );
313 
314  kEqn.ref().relax();
315  fvOptions.constrain(kEqn.ref());
316  solve(kEqn);
317  fvOptions.correct(k_);
318  bound(k_, this->kMin_);
319 
320  correctNut();
321 }
322 
323 
324 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
325 
326 } // End namespace RASModels
327 } // End namespace Foam
328 
329 // ************************************************************************* //
surfaceScalarField & phi
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
U
Definition: pEqn.H:83
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:51
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:241
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:100
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite-volume options.
Definition: fvOptions.H:52
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:101
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
fv::options & fvOptions
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:66
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
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
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-5.0/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
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:485
Renormalization group k-epsilon turbulence model for incompressible and compressible flows...
Definition: RNGkEpsilon.H:85
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:218
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
#define R(A, B, C, D, E, F, K, M)
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
dimensioned< scalar > mag(const dimensioned< Type > &)
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:103
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:92
scalar nut
Namespace for OpenFOAM.
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:99