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-2023 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  tmp<volScalarField> tCmuk2(Cmu_*sqr(k_));
44  epsilon_ = max(epsilon_, tCmuk2()/(this->nutMaxCoeff_*this->nu()));
45  return tCmuk2;
46 }
47 
48 
49 template<class BasicMomentumTransportModel>
51 {
52  this->nut_ = boundEpsilon()/epsilon_;
53  this->nut_.correctBoundaryConditions();
54  fvConstraints::New(this->mesh_).constrain(this->nut_);
55 }
56 
57 
58 template<class BasicMomentumTransportModel>
60 {
61  return tmp<fvScalarMatrix>
62  (
63  new fvScalarMatrix
64  (
65  k_,
66  dimVolume*this->rho_.dimensions()*k_.dimensions()
67  /dimTime
68  )
69  );
70 }
71 
72 
73 template<class BasicMomentumTransportModel>
76 {
77  return tmp<fvScalarMatrix>
78  (
79  new fvScalarMatrix
80  (
81  epsilon_,
82  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
83  /dimTime
84  )
85  );
86 }
87 
88 
89 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
90 
91 template<class BasicMomentumTransportModel>
93 (
94  const alphaField& alpha,
95  const rhoField& rho,
96  const volVectorField& U,
97  const surfaceScalarField& alphaRhoPhi,
98  const surfaceScalarField& phi,
99  const viscosity& viscosity,
100  const word& type
101 )
102 :
103  eddyViscosity<RASModel<BasicMomentumTransportModel>>
104  (
105  type,
106  alpha,
107  rho,
108  U,
109  alphaRhoPhi,
110  phi,
111  viscosity
112  ),
113 
114  Cmu_
115  (
116  dimensioned<scalar>::lookupOrAddToDict
117  (
118  "Cmu",
119  this->coeffDict_,
120  0.0845
121  )
122  ),
123  C1_
124  (
125  dimensioned<scalar>::lookupOrAddToDict
126  (
127  "C1",
128  this->coeffDict_,
129  1.42
130  )
131  ),
132  C2_
133  (
134  dimensioned<scalar>::lookupOrAddToDict
135  (
136  "C2",
137  this->coeffDict_,
138  1.68
139  )
140  ),
141  C3_
142  (
143  dimensioned<scalar>::lookupOrAddToDict
144  (
145  "C3",
146  this->coeffDict_,
147  0
148  )
149  ),
150  sigmak_
151  (
152  dimensioned<scalar>::lookupOrAddToDict
153  (
154  "sigmak",
155  this->coeffDict_,
156  0.71942
157  )
158  ),
159  sigmaEps_
160  (
161  dimensioned<scalar>::lookupOrAddToDict
162  (
163  "sigmaEps",
164  this->coeffDict_,
165  0.71942
166  )
167  ),
168  eta0_
169  (
170  dimensioned<scalar>::lookupOrAddToDict
171  (
172  "eta0",
173  this->coeffDict_,
174  4.38
175  )
176  ),
177  beta_
178  (
179  dimensioned<scalar>::lookupOrAddToDict
180  (
181  "beta",
182  this->coeffDict_,
183  0.012
184  )
185  ),
186 
187  k_
188  (
189  IOobject
190  (
191  this->groupName("k"),
192  this->runTime_.name(),
193  this->mesh_,
194  IOobject::MUST_READ,
195  IOobject::AUTO_WRITE
196  ),
197  this->mesh_
198  ),
199  epsilon_
200  (
201  IOobject
202  (
203  this->groupName("epsilon"),
204  this->runTime_.name(),
205  this->mesh_,
206  IOobject::MUST_READ,
207  IOobject::AUTO_WRITE
208  ),
209  this->mesh_
210  )
211 {
212  bound(k_, this->kMin_);
213  boundEpsilon();
214 
215  if (type == typeName)
216  {
217  this->printCoeffs(type);
218  }
219 }
220 
221 
222 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
223 
224 template<class BasicMomentumTransportModel>
226 {
228  {
229  Cmu_.readIfPresent(this->coeffDict());
230  C1_.readIfPresent(this->coeffDict());
231  C2_.readIfPresent(this->coeffDict());
232  C3_.readIfPresent(this->coeffDict());
233  sigmak_.readIfPresent(this->coeffDict());
234  sigmaEps_.readIfPresent(this->coeffDict());
235  eta0_.readIfPresent(this->coeffDict());
236  beta_.readIfPresent(this->coeffDict());
237 
238  return true;
239  }
240  else
241  {
242  return false;
243  }
244 }
245 
246 
247 template<class BasicMomentumTransportModel>
249 {
250  if (!this->turbulence_)
251  {
252  return;
253  }
254 
255  // Local references
256  const alphaField& alpha = this->alpha_;
257  const rhoField& rho = this->rho_;
258  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
259  const volVectorField& U = this->U_;
260  volScalarField& nut = this->nut_;
261  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
263  (
264  Foam::fvConstraints::New(this->mesh_)
265  );
266 
268 
270  (
271  typedName("divU"),
272  fvc::div(fvc::absolute(this->phi(), U))()
273  );
274 
275  tmp<volTensorField> tgradU = fvc::grad(U);
277  (
278  typedName("S2"),
279  (tgradU().v() && dev(twoSymm(tgradU().v())))
280  );
281  tgradU.clear();
282 
283  volScalarField::Internal G(this->GName(), nut()*S2);
284 
286  (
287  typedName("eta"),
288  sqrt(mag(S2))*k_()/epsilon_()
289  );
290 
291  volScalarField::Internal eta3(typedName("eta3"), eta*sqr(eta));
292 
294  (
295  typedName("R"),
296  ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
297  );
298 
299  // Update epsilon and G at the wall
300  epsilon_.boundaryFieldRef().updateCoeffs();
301 
302  // Dissipation equation
303  tmp<fvScalarMatrix> epsEqn
304  (
305  fvm::ddt(alpha, rho, epsilon_)
306  + fvm::div(alphaRhoPhi, epsilon_)
307  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
308  ==
309  (C1_ - R)*alpha()*rho()*G*epsilon_()/k_()
310  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
311  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
312  + epsilonSource()
313  + fvModels.source(alpha, rho, epsilon_)
314  );
315 
316  epsEqn.ref().relax();
317  fvConstraints.constrain(epsEqn.ref());
318  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
319  solve(epsEqn);
320  fvConstraints.constrain(epsilon_);
321  boundEpsilon();
322 
323 
324  // Turbulent kinetic energy equation
325 
327  (
328  fvm::ddt(alpha, rho, k_)
329  + fvm::div(alphaRhoPhi, k_)
330  - fvm::laplacian(alpha*rho*DkEff(), k_)
331  ==
332  alpha()*rho()*G
333  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
334  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
335  + kSource()
336  + fvModels.source(alpha, rho, k_)
337  );
338 
339  kEqn.ref().relax();
340  fvConstraints.constrain(kEqn.ref());
341  solve(kEqn);
343  bound(k_, this->kMin_);
344 
345  correctNut();
346 }
347 
348 
349 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
350 
351 } // End namespace RASModels
352 } // End namespace Foam
353 
354 // ************************************************************************* //
Bound the given scalar field where it is below the specified minimum.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
virtual tmp< fvScalarMatrix > epsilonSource() const
Source term for the epsilon equation.
Definition: RNGkEpsilon.C:75
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: RNGkEpsilon.C:248
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
Definition: RNGkEpsilon.C:41
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: RNGkEpsilon.C:50
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
Definition: RNGkEpsilon.C:59
RNGkEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
Definition: RNGkEpsilon.C:93
virtual bool read()
Re-read model coefficients if they have changed.
Definition: RNGkEpsilon.C:225
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite volume constraints.
Definition: fvConstraints.H:67
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Finite volume models.
Definition: fvModels.H:65
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
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
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
const scalar nut
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor sqr(const dimensionedVector &dv)
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.