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-2024 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_("Cmu", this->coeffDict(), 0.0845),
115  C1_("C1", this->coeffDict(), 1.42),
116  C2_("C2", this->coeffDict(), 1.68),
117  C3_("C3", this->coeffDict(), 0),
118  sigmak_("sigmak", this->coeffDict(), 0.71942),
119  sigmaEps_("sigmaEps", this->coeffDict(), 0.71942),
120  eta0_("eta0", this->coeffDict(), 4.38),
121  beta_("beta", this->coeffDict(), 0.012),
122 
123  k_
124  (
125  IOobject
126  (
127  this->groupName("k"),
128  this->runTime_.name(),
129  this->mesh_,
130  IOobject::MUST_READ,
131  IOobject::AUTO_WRITE
132  ),
133  this->mesh_
134  ),
135  epsilon_
136  (
137  IOobject
138  (
139  this->groupName("epsilon"),
140  this->runTime_.name(),
141  this->mesh_,
142  IOobject::MUST_READ,
143  IOobject::AUTO_WRITE
144  ),
145  this->mesh_
146  )
147 {
148  bound(k_, this->kMin_);
149  boundEpsilon();
150 }
151 
152 
153 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
154 
155 template<class BasicMomentumTransportModel>
157 {
159  {
160  Cmu_.readIfPresent(this->coeffDict());
161  C1_.readIfPresent(this->coeffDict());
162  C2_.readIfPresent(this->coeffDict());
163  C3_.readIfPresent(this->coeffDict());
164  sigmak_.readIfPresent(this->coeffDict());
165  sigmaEps_.readIfPresent(this->coeffDict());
166  eta0_.readIfPresent(this->coeffDict());
167  beta_.readIfPresent(this->coeffDict());
168 
169  return true;
170  }
171  else
172  {
173  return false;
174  }
175 }
176 
177 
178 template<class BasicMomentumTransportModel>
180 {
181  if (!this->turbulence_)
182  {
183  return;
184  }
185 
186  // Local references
187  const alphaField& alpha = this->alpha_;
188  const rhoField& rho = this->rho_;
189  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
190  const volVectorField& U = this->U_;
191  volScalarField& nut = this->nut_;
192  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
194  (
195  Foam::fvConstraints::New(this->mesh_)
196  );
197 
199 
201  (
202  typedName("divU"),
203  fvc::div(fvc::absolute(this->phi(), U))()
204  );
205 
206  tmp<volTensorField> tgradU = fvc::grad(U);
208  (
209  typedName("S2"),
210  (tgradU().v() && dev(twoSymm(tgradU().v())))
211  );
212  tgradU.clear();
213 
214  volScalarField::Internal G(this->GName(), nut()*S2);
215 
217  (
218  typedName("eta"),
219  sqrt(mag(S2))*k_()/epsilon_()
220  );
221 
222  volScalarField::Internal eta3(typedName("eta3"), eta*sqr(eta));
223 
225  (
226  typedName("R"),
227  ((eta*(-eta/eta0_ + scalar(1)))/(beta_*eta3 + scalar(1)))
228  );
229 
230  // Update epsilon and G at the wall
231  epsilon_.boundaryFieldRef().updateCoeffs();
232 
233  // Dissipation equation
234  tmp<fvScalarMatrix> epsEqn
235  (
236  fvm::ddt(alpha, rho, epsilon_)
237  + fvm::div(alphaRhoPhi, epsilon_)
238  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
239  ==
240  (C1_ - R)*alpha()*rho()*G*epsilon_()/k_()
241  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
242  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
243  + epsilonSource()
244  + fvModels.source(alpha, rho, epsilon_)
245  );
246 
247  epsEqn.ref().relax();
248  fvConstraints.constrain(epsEqn.ref());
249  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
250  solve(epsEqn);
251  fvConstraints.constrain(epsilon_);
252  boundEpsilon();
253 
254 
255  // Turbulent kinetic energy equation
256 
258  (
259  fvm::ddt(alpha, rho, k_)
260  + fvm::div(alphaRhoPhi, k_)
261  - fvm::laplacian(alpha*rho*DkEff(), k_)
262  ==
263  alpha()*rho()*G
264  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
265  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
266  + kSource()
267  + fvModels.source(alpha, rho, k_)
268  );
269 
270  kEqn.ref().relax();
271  fvConstraints.constrain(kEqn.ref());
272  solve(kEqn);
274  bound(k_, this->kMin_);
275 
276  correctNut();
277 }
278 
279 
280 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
281 
282 } // End namespace RASModels
283 } // End namespace Foam
284 
285 // ************************************************************************* //
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:179
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:156
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:103
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:253
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
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 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.
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
const dimensionSet dimVolume
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:181
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
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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.