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-2020 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 {
56  return tmp<fvScalarMatrix>
57  (
58  new fvScalarMatrix
59  (
60  k_,
61  dimVolume*this->rho_.dimensions()*k_.dimensions()
62  /dimTime
63  )
64  );
65 }
66 
67 
68 template<class BasicTurbulenceModel>
70 {
71  return tmp<fvScalarMatrix>
72  (
73  new fvScalarMatrix
74  (
75  epsilon_,
76  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
77  /dimTime
78  )
79  );
80 }
81 
82 
83 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
84 
85 template<class BasicTurbulenceModel>
87 (
88  const alphaField& alpha,
89  const rhoField& rho,
90  const volVectorField& U,
91  const surfaceScalarField& alphaRhoPhi,
92  const surfaceScalarField& phi,
93  const transportModel& transport,
94  const word& propertiesName,
95  const word& type
96 )
97 :
99  (
100  type,
101  alpha,
102  rho,
103  U,
104  alphaRhoPhi,
105  phi,
106  transport,
107  propertiesName
108  ),
109 
110  Cmu_
111  (
113  (
114  "Cmu",
115  this->coeffDict_,
116  0.0845
117  )
118  ),
119  C1_
120  (
122  (
123  "C1",
124  this->coeffDict_,
125  1.42
126  )
127  ),
128  C2_
129  (
131  (
132  "C2",
133  this->coeffDict_,
134  1.68
135  )
136  ),
137  C3_
138  (
140  (
141  "C3",
142  this->coeffDict_,
143  0
144  )
145  ),
146  sigmak_
147  (
149  (
150  "sigmak",
151  this->coeffDict_,
152  0.71942
153  )
154  ),
155  sigmaEps_
156  (
158  (
159  "sigmaEps",
160  this->coeffDict_,
161  0.71942
162  )
163  ),
164  eta0_
165  (
167  (
168  "eta0",
169  this->coeffDict_,
170  4.38
171  )
172  ),
173  beta_
174  (
176  (
177  "beta",
178  this->coeffDict_,
179  0.012
180  )
181  ),
182 
183  k_
184  (
185  IOobject
186  (
187  IOobject::groupName("k", alphaRhoPhi.group()),
188  this->runTime_.timeName(),
189  this->mesh_,
192  ),
193  this->mesh_
194  ),
195  epsilon_
196  (
197  IOobject
198  (
199  IOobject::groupName("epsilon", alphaRhoPhi.group()),
200  this->runTime_.timeName(),
201  this->mesh_,
204  ),
205  this->mesh_
206  )
207 {
208  bound(k_, this->kMin_);
209  bound(epsilon_, this->epsilonMin_);
210 
211  if (type == typeName)
212  {
213  this->printCoeffs(type);
214  }
215 }
216 
217 
218 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
219 
220 template<class BasicTurbulenceModel>
222 {
224  {
225  Cmu_.readIfPresent(this->coeffDict());
226  C1_.readIfPresent(this->coeffDict());
227  C2_.readIfPresent(this->coeffDict());
228  C3_.readIfPresent(this->coeffDict());
229  sigmak_.readIfPresent(this->coeffDict());
230  sigmaEps_.readIfPresent(this->coeffDict());
231  eta0_.readIfPresent(this->coeffDict());
232  beta_.readIfPresent(this->coeffDict());
233 
234  return true;
235  }
236 
237  return false;
238 }
239 
240 
241 template<class BasicTurbulenceModel>
243 {
244  if (!this->turbulence_)
245  {
246  return;
247  }
248 
249  // Local references
250  const alphaField& alpha = this->alpha_;
251  const rhoField& rho = this->rho_;
252  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
253  const volVectorField& U = this->U_;
254  const volScalarField& nut = this->nut_;
255  fv::options& fvOptions(fv::options::New(this->mesh_));
256 
258 
260  (
261  fvc::div(fvc::absolute(this->phi(), U))().v()
262  );
263 
264  tmp<volTensorField> tgradU = fvc::grad(U);
265  const volScalarField::Internal GbyNu
266  (
267  tgradU().v() && dev(twoSymm(tgradU().v()))
268  );
269  tgradU.clear();
270 
271  const volScalarField::Internal G(this->GName(), nut()*GbyNu);
272 
273  const volScalarField::Internal eta(sqrt(mag(GbyNu))*k_/epsilon_);
274  const volScalarField::Internal eta3(eta*sqr(eta));
275 
277  (
278  ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
279  );
280 
281  // Update epsilon and G at the wall
282  epsilon_.boundaryFieldRef().updateCoeffs();
283 
284  // Dissipation equation
285  tmp<fvScalarMatrix> epsEqn
286  (
287  fvm::ddt(alpha, rho, epsilon_)
288  + fvm::div(alphaRhoPhi, epsilon_)
289  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
290  ==
291  (C1_ - R)*alpha()*rho()*GbyNu*Cmu_*k_()
292  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
293  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
294  + epsilonSource()
295  + fvOptions(alpha, rho, epsilon_)
296  );
297 
298  epsEqn.ref().relax();
299  fvOptions.constrain(epsEqn.ref());
300  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
301  solve(epsEqn);
302  fvOptions.correct(epsilon_);
303  bound(epsilon_, this->epsilonMin_);
304 
305 
306  // Turbulent kinetic energy equation
307 
309  (
310  fvm::ddt(alpha, rho, k_)
311  + fvm::div(alphaRhoPhi, k_)
312  - fvm::laplacian(alpha*rho*DkEff(), k_)
313  ==
314  alpha()*rho()*G
315  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
316  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
317  + kSource()
318  + fvOptions(alpha, rho, k_)
319  );
320 
321  kEqn.ref().relax();
322  fvOptions.constrain(kEqn.ref());
323  solve(kEqn);
324  fvOptions.correct(k_);
325  bound(k_, this->kMin_);
326 
327  correctNut();
328 }
329 
330 
331 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
332 
333 } // End namespace RASModels
334 } // End namespace Foam
335 
336 // ************************************************************************* //
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
Foam::IOobject::AUTO_WRITE
Definition: IOobject.H:194
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:54
fvOptions.H
Foam::word
A class for handling words, derived from Foam::string.
Definition: word.H:65
Foam::fv::optionList::correct
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
Definition: fvOptionListTemplates.C:355
Foam::tmp::clear
void clear() const noexcept
Definition: tmpI.H:287
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:52
Foam::fv::options::New
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:103
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
Foam::RASModels::RNGkEpsilon::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: RNGkEpsilon.C:54
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:35
Foam::read
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:108
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:49
Foam::RASModels::RNGkEpsilon::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: RNGkEpsilon.C:69
rho
rho
Definition: readInitialConditions.H:88
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::RASModels::RNGkEpsilon
Renormalization group k-epsilon turbulence model for incompressible and compressible flows.
Definition: RNGkEpsilon.H:88
R
#define R(A, B, C, D, E, F, K, M)
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::RASModels::RNGkEpsilon::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:242
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
fvOptions
fv::options & fvOptions
Definition: setRegionFluidFields.H:23
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:48
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::solve
SolverPerformance< Type > solve(faMatrix< Type > &, Istream &)
Solve returning the solution statistics given convergence tolerance.
divU
zeroField divU
Definition: alphaSuSp.H:3
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::fv::options
Finite-volume options.
Definition: fvOptions.H:55
Foam::RASModel::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: RASModel.H:98
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:42
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::RASModel::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: RASModel.H:97
Foam::RASModels::RNGkEpsilon::correctNut
virtual void correctNut()
Definition: RNGkEpsilon.C:43
Foam::RASModel::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: RASModel.H:99
U
U
Definition: pEqn.H:72
Foam::RASModels::RNGkEpsilon::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:221
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:48
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::mag
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
Foam::type
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: MSwindows.C:590
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:55
RNGkEpsilon.H
Foam::IOobject::groupName
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
nut
scalar nut
Definition: evaluateNearWall.H:5
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::GeometricField< vector, fvPatchField, volMesh >
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:48
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:95
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:190
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
Foam::IOobject::MUST_READ
Definition: IOobject.H:185
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:106