RNGkEpsilon.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2011-2020 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 BasicMomentumTransportModel>
41 {
42  this->nut_ = Cmu_*sqr(k_)/epsilon_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 }
46 
47 
48 template<class BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
66 {
67  return tmp<fvScalarMatrix>
68  (
69  new fvScalarMatrix
70  (
71  epsilon_,
72  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
73  /dimTime
74  )
75  );
76 }
77 
78 
79 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
80 
81 template<class BasicMomentumTransportModel>
83 (
84  const alphaField& alpha,
85  const rhoField& rho,
86  const volVectorField& U,
87  const surfaceScalarField& alphaRhoPhi,
88  const surfaceScalarField& phi,
89  const transportModel& transport,
90  const word& type
91 )
92 :
94  (
95  type,
96  alpha,
97  rho,
98  U,
99  alphaRhoPhi,
100  phi,
101  transport
102  ),
103 
104  Cmu_
105  (
107  (
108  "Cmu",
109  this->coeffDict_,
110  0.0845
111  )
112  ),
113  C1_
114  (
116  (
117  "C1",
118  this->coeffDict_,
119  1.42
120  )
121  ),
122  C2_
123  (
125  (
126  "C2",
127  this->coeffDict_,
128  1.68
129  )
130  ),
131  C3_
132  (
134  (
135  "C3",
136  this->coeffDict_,
137  0
138  )
139  ),
140  sigmak_
141  (
143  (
144  "sigmak",
145  this->coeffDict_,
146  0.71942
147  )
148  ),
149  sigmaEps_
150  (
152  (
153  "sigmaEps",
154  this->coeffDict_,
155  0.71942
156  )
157  ),
158  eta0_
159  (
161  (
162  "eta0",
163  this->coeffDict_,
164  4.38
165  )
166  ),
167  beta_
168  (
170  (
171  "beta",
172  this->coeffDict_,
173  0.012
174  )
175  ),
176 
177  k_
178  (
179  IOobject
180  (
181  IOobject::groupName("k", alphaRhoPhi.group()),
182  this->runTime_.timeName(),
183  this->mesh_,
186  ),
187  this->mesh_
188  ),
189  epsilon_
190  (
191  IOobject
192  (
193  IOobject::groupName("epsilon", alphaRhoPhi.group()),
194  this->runTime_.timeName(),
195  this->mesh_,
198  ),
199  this->mesh_
200  )
201 {
202  bound(k_, this->kMin_);
203  bound(epsilon_, this->epsilonMin_);
204 
205  if (type == typeName)
206  {
207  this->printCoeffs(type);
208  }
209 }
210 
211 
212 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
213 
214 template<class BasicMomentumTransportModel>
216 {
218  {
219  Cmu_.readIfPresent(this->coeffDict());
220  C1_.readIfPresent(this->coeffDict());
221  C2_.readIfPresent(this->coeffDict());
222  C3_.readIfPresent(this->coeffDict());
223  sigmak_.readIfPresent(this->coeffDict());
224  sigmaEps_.readIfPresent(this->coeffDict());
225  eta0_.readIfPresent(this->coeffDict());
226  beta_.readIfPresent(this->coeffDict());
227 
228  return true;
229  }
230  else
231  {
232  return false;
233  }
234 }
235 
236 
237 template<class BasicMomentumTransportModel>
239 {
240  if (!this->turbulence_)
241  {
242  return;
243  }
244 
245  // Local references
246  const alphaField& alpha = this->alpha_;
247  const rhoField& rho = this->rho_;
248  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
249  const volVectorField& U = this->U_;
250  volScalarField& nut = this->nut_;
251  fv::options& fvOptions(fv::options::New(this->mesh_));
252 
254 
256  (
257  modelName("divU"),
258  fvc::div(fvc::absolute(this->phi(), U))()
259  );
260 
261  tmp<volTensorField> tgradU = fvc::grad(U);
263  (
264  modelName("S2"),
265  (tgradU().v() && dev(twoSymm(tgradU().v())))
266  );
267  tgradU.clear();
268 
269  volScalarField::Internal G(this->GName(), nut()*S2);
270 
272  (
273  modelName("eta"),
274  sqrt(mag(S2))*k_()/epsilon_()
275  );
276 
277  volScalarField::Internal eta3(modelName("eta3"), eta*sqr(eta));
278 
280  (
281  modelName("R"),
282  ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
283  );
284 
285  // Update epsilon and G at the wall
286  epsilon_.boundaryFieldRef().updateCoeffs();
287 
288  // Dissipation equation
289  tmp<fvScalarMatrix> epsEqn
290  (
291  fvm::ddt(alpha, rho, epsilon_)
292  + fvm::div(alphaRhoPhi, epsilon_)
293  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
294  ==
295  (C1_ - R)*alpha()*rho()*G*epsilon_()/k_()
296  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
297  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
298  + epsilonSource()
299  + fvOptions(alpha, rho, epsilon_)
300  );
301 
302  epsEqn.ref().relax();
303  fvOptions.constrain(epsEqn.ref());
304  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
305  solve(epsEqn);
306  fvOptions.correct(epsilon_);
307  bound(epsilon_, this->epsilonMin_);
308 
309 
310  // Turbulent kinetic energy equation
311 
313  (
314  fvm::ddt(alpha, rho, k_)
315  + fvm::div(alphaRhoPhi, k_)
316  - fvm::laplacian(alpha*rho*DkEff(), k_)
317  ==
318  alpha()*rho()*G
319  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
320  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
321  + kSource()
322  + fvOptions(alpha, rho, k_)
323  );
324 
325  kEqn.ref().relax();
326  fvOptions.constrain(kEqn.ref());
327  solve(kEqn);
328  fvOptions.correct(k_);
329  bound(k_, this->kMin_);
330 
331  correctNut();
332 }
333 
334 
335 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
336 
337 } // End namespace RASModels
338 } // End namespace Foam
339 
340 // ************************************************************************* //
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
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
fv::options & fvOptions
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
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:49
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:238
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:49
Finite-volume options.
Definition: fvOptions.H:52
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:88
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
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:65
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
RNGkEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: RNGkEpsilon.C:83
BasicMomentumTransportModel::transportModel transportModel
Definition: RASModel.H:90
A class for handling words, derived from string.
Definition: word.H:59
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(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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:68
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:215
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:89
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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:101
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
scalar nut
Namespace for OpenFOAM.