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-2021 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 "fvModels.H"
28 #include "fvConstraints.H"
29 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
42 {
43  this->nut_ = Cmu_*sqr(k_)/epsilon_;
44  this->nut_.correctBoundaryConditions();
45  fvConstraints::New(this->mesh_).constrain(this->nut_);
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  return tmp<fvScalarMatrix>
53  (
54  new fvScalarMatrix
55  (
56  k_,
57  dimVolume*this->rho_.dimensions()*k_.dimensions()
58  /dimTime
59  )
60  );
61 }
62 
63 
64 template<class BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
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& type
92 )
93 :
95  (
96  type,
97  alpha,
98  rho,
99  U,
100  alphaRhoPhi,
101  phi,
102  transport
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
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", alphaRhoPhi.group()),
183  this->runTime_.timeName(),
184  this->mesh_,
187  ),
188  this->mesh_
189  ),
190  epsilon_
191  (
192  IOobject
193  (
194  IOobject::groupName("epsilon", alphaRhoPhi.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 BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
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  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
254  (
255  Foam::fvConstraints::New(this->mesh_)
256  );
257 
259 
261  (
262  modelName("divU"),
263  fvc::div(fvc::absolute(this->phi(), U))()
264  );
265 
266  tmp<volTensorField> tgradU = fvc::grad(U);
268  (
269  modelName("S2"),
270  (tgradU().v() && dev(twoSymm(tgradU().v())))
271  );
272  tgradU.clear();
273 
274  volScalarField::Internal G(this->GName(), nut()*S2);
275 
277  (
278  modelName("eta"),
279  sqrt(mag(S2))*k_()/epsilon_()
280  );
281 
282  volScalarField::Internal eta3(modelName("eta3"), eta*sqr(eta));
283 
285  (
286  modelName("R"),
287  ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
288  );
289 
290  // Update epsilon and G at the wall
291  epsilon_.boundaryFieldRef().updateCoeffs();
292 
293  // Dissipation equation
294  tmp<fvScalarMatrix> epsEqn
295  (
296  fvm::ddt(alpha, rho, epsilon_)
297  + fvm::div(alphaRhoPhi, epsilon_)
298  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
299  ==
300  (C1_ - R)*alpha()*rho()*G*epsilon_()/k_()
301  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
302  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
303  + epsilonSource()
304  + fvModels.source(alpha, rho, epsilon_)
305  );
306 
307  epsEqn.ref().relax();
308  fvConstraints.constrain(epsEqn.ref());
309  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
310  solve(epsEqn);
311  fvConstraints.constrain(epsilon_);
312  bound(epsilon_, this->epsilonMin_);
313 
314 
315  // Turbulent kinetic energy equation
316 
318  (
319  fvm::ddt(alpha, rho, k_)
320  + fvm::div(alphaRhoPhi, k_)
321  - fvm::laplacian(alpha*rho*DkEff(), k_)
322  ==
323  alpha()*rho()*G
324  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
325  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
326  + kSource()
327  + fvModels.source(alpha, rho, k_)
328  );
329 
330  kEqn.ref().relax();
331  fvConstraints.constrain(kEqn.ref());
332  solve(kEqn);
333  fvConstraints.constrain(k_);
334  bound(k_, this->kMin_);
335 
336  correctNut();
337 }
338 
339 
340 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
341 
342 } // End namespace RASModels
343 } // End namespace Foam
344 
345 // ************************************************************************* //
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:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
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:50
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:239
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:88
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
const dimensionSet dimTime
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:66
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:84
BasicMomentumTransportModel::transportModel transportModel
Definition: RASModel.H:90
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/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
static word groupName(Name name, const word &group)
Finite volume models.
Definition: fvModels.H:60
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
phi
Definition: correctPhi.H:3
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.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:216
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
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)
Finite volume constraints.
Definition: fvConstraints.H:57
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 dimVolume
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
dimensioned< scalar > mag(const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
scalar nut
Namespace for OpenFOAM.