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 | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "RNGkEpsilon.H"
30 #include "fvOptions.H"
31 #include "bound.H"
32 
33 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace RASModels
38 {
39 
40 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
41 
42 template<class BasicTurbulenceModel>
44 {
45  this->nut_ = Cmu_*sqr(k_)/epsilon_;
46  this->nut_.correctBoundaryConditions();
47  fv::options::New(this->mesh_).correct(this->nut_);
48 
49  BasicTurbulenceModel::correctNut();
50 }
51 
52 
53 template<class BasicTurbulenceModel>
55 {
57  (
58  k_,
59  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
60  );
61 }
62 
63 
64 template<class BasicTurbulenceModel>
66 {
68  (
69  epsilon_,
70  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()/dimTime
71  );
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
77 template<class BasicTurbulenceModel>
79 (
80  const alphaField& alpha,
81  const rhoField& rho,
82  const volVectorField& U,
83  const surfaceScalarField& alphaRhoPhi,
84  const surfaceScalarField& phi,
85  const transportModel& transport,
86  const word& propertiesName,
87  const word& type
88 )
89 :
90  eddyViscosity<RASModel<BasicTurbulenceModel>>
91  (
92  type,
93  alpha,
94  rho,
95  U,
96  alphaRhoPhi,
97  phi,
98  transport,
99  propertiesName
100  ),
101 
102  Cmu_
103  (
104  dimensioned<scalar>::getOrAddToDict
105  (
106  "Cmu",
107  this->coeffDict_,
108  0.0845
109  )
110  ),
111  C1_
112  (
113  dimensioned<scalar>::getOrAddToDict
114  (
115  "C1",
116  this->coeffDict_,
117  1.42
118  )
119  ),
120  C2_
121  (
122  dimensioned<scalar>::getOrAddToDict
123  (
124  "C2",
125  this->coeffDict_,
126  1.68
127  )
128  ),
129  C3_
130  (
131  dimensioned<scalar>::getOrAddToDict
132  (
133  "C3",
134  this->coeffDict_,
135  0
136  )
137  ),
138  sigmak_
139  (
140  dimensioned<scalar>::getOrAddToDict
141  (
142  "sigmak",
143  this->coeffDict_,
144  0.71942
145  )
146  ),
147  sigmaEps_
148  (
149  dimensioned<scalar>::getOrAddToDict
150  (
151  "sigmaEps",
152  this->coeffDict_,
153  0.71942
154  )
155  ),
156  eta0_
157  (
158  dimensioned<scalar>::getOrAddToDict
159  (
160  "eta0",
161  this->coeffDict_,
162  4.38
163  )
164  ),
165  beta_
166  (
167  dimensioned<scalar>::getOrAddToDict
168  (
169  "beta",
170  this->coeffDict_,
171  0.012
172  )
173  ),
174 
175  k_
176  (
177  IOobject
178  (
179  IOobject::groupName("k", alphaRhoPhi.group()),
180  this->runTime_.timeName(),
181  this->mesh_,
182  IOobject::MUST_READ,
183  IOobject::AUTO_WRITE
184  ),
185  this->mesh_
186  ),
187  epsilon_
188  (
189  IOobject
190  (
191  IOobject::groupName("epsilon", alphaRhoPhi.group()),
192  this->runTime_.timeName(),
193  this->mesh_,
194  IOobject::MUST_READ,
195  IOobject::AUTO_WRITE
196  ),
197  this->mesh_
198  )
199 {
200  bound(k_, this->kMin_);
201  bound(epsilon_, this->epsilonMin_);
202 
203  if (type == typeName)
204  {
205  this->printCoeffs(type);
206  }
207 }
208 
209 
210 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
211 
212 template<class BasicTurbulenceModel>
214 {
216  {
217  Cmu_.readIfPresent(this->coeffDict());
218  C1_.readIfPresent(this->coeffDict());
219  C2_.readIfPresent(this->coeffDict());
220  C3_.readIfPresent(this->coeffDict());
221  sigmak_.readIfPresent(this->coeffDict());
222  sigmaEps_.readIfPresent(this->coeffDict());
223  eta0_.readIfPresent(this->coeffDict());
224  beta_.readIfPresent(this->coeffDict());
225 
226  return true;
227  }
228 
229  return false;
230 }
231 
232 
233 template<class BasicTurbulenceModel>
235 {
236  if (!this->turbulence_)
237  {
238  return;
239  }
240 
241  // Local references
242  const alphaField& alpha = this->alpha_;
243  const rhoField& rho = this->rho_;
244  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
245  const volVectorField& U = this->U_;
246  const volScalarField& nut = this->nut_;
247  fv::options& fvOptions(fv::options::New(this->mesh_));
248 
249  eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
250 
252  (
253  fvc::div(fvc::absolute(this->phi(), U))().v()
254  );
255 
256  tmp<volTensorField> tgradU = fvc::grad(U);
257  const volScalarField::Internal GbyNu
258  (
259  IOobject::scopedName(this->type(), "GbyNu"),
260  tgradU().v() && devTwoSymm(tgradU().v())
261  );
262  tgradU.clear();
263 
264  const volScalarField::Internal G(this->GName(), nut()*GbyNu);
265 
266  const volScalarField::Internal eta(sqrt(mag(GbyNu))*k_/epsilon_);
267  const volScalarField::Internal 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  // Push any changed cell values to coupled neighbours
277  epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
278 
279  // Dissipation equation
280  tmp<fvScalarMatrix> epsEqn
281  (
282  fvm::ddt(alpha, rho, epsilon_)
283  + fvm::div(alphaRhoPhi, epsilon_)
284  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
285  ==
286  (C1_ - R)*alpha()*rho()*GbyNu*Cmu_*k_()
287  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
288  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
289  + epsilonSource()
290  + fvOptions(alpha, rho, epsilon_)
291  );
292 
293  epsEqn.ref().relax();
294  fvOptions.constrain(epsEqn.ref());
295  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
296  solve(epsEqn);
297  fvOptions.correct(epsilon_);
298  bound(epsilon_, this->epsilonMin_);
299 
300 
301  // Turbulent kinetic energy equation
302 
303  tmp<fvScalarMatrix> kEqn
304  (
305  fvm::ddt(alpha, rho, k_)
306  + fvm::div(alphaRhoPhi, k_)
307  - fvm::laplacian(alpha*rho*DkEff(), k_)
308  ==
309  alpha()*rho()*G
310  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
311  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
312  + kSource()
313  + fvOptions(alpha, rho, k_)
314  );
315 
316  kEqn.ref().relax();
317  fvOptions.constrain(kEqn.ref());
318  solve(kEqn);
319  fvOptions.correct(k_);
320  bound(k_, this->kMin_);
321 
322  correctNut();
323 }
324 
325 
326 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
327 
328 } // End namespace RASModels
329 } // End namespace Foam
330 
331 // ************************************************************************* //
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:109
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
const dimensionedScalar G
Newtonian constant of gravitation.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:47
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:227
dimensionedScalar sqrt(const dimensionedScalar &ds)
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
DimensionedField< scalar, volMesh > Internal
The internal field type from which this GeometricField is derived.
fv::options & fvOptions
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:77
constexpr const char *const group
Group name for atomic constants.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:72
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:58
static word scopedName(const std::string &scope, const word &name)
Create scope:name or scope_name string.
Definition: IOobjectI.H:40
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:27
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
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/chef2/andy/OpenFOAM/release/v2406/OpenFOAM-v2406/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=clamp((rho - rholSat)/(rhovSat - rholSat), zero_one{});alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:63
Renormalization group k-epsilon turbulence model for incompressible and compressible flows...
Definition: RNGkEpsilon.H:84
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:206
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:183
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:29
#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:41
Base-class for all transport models used by the incompressible turbulence models. ...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:180
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:114
scalar nut
Namespace for OpenFOAM.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491