kEpsilon.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 "kEpsilon.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>
75 {
76  return tmp<fvScalarMatrix>
77  (
78  new fvScalarMatrix
79  (
80  epsilon_,
81  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
82  /dimTime
83  )
84  );
85 }
86 
87 
88 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
89 
90 template<class BasicMomentumTransportModel>
92 (
93  const alphaField& alpha,
94  const rhoField& rho,
95  const volVectorField& U,
96  const surfaceScalarField& alphaRhoPhi,
97  const surfaceScalarField& phi,
98  const viscosity& viscosity,
99  const word& type
100 )
101 :
102  eddyViscosity<RASModel<BasicMomentumTransportModel>>
103  (
104  type,
105  alpha,
106  rho,
107  U,
108  alphaRhoPhi,
109  phi,
110  viscosity
111  ),
112 
113  Cmu_
114  (
115  dimensioned<scalar>::lookupOrAddToDict
116  (
117  "Cmu",
118  this->coeffDict_,
119  0.09
120  )
121  ),
122  C1_
123  (
124  dimensioned<scalar>::lookupOrAddToDict
125  (
126  "C1",
127  this->coeffDict_,
128  1.44
129  )
130  ),
131  C2_
132  (
133  dimensioned<scalar>::lookupOrAddToDict
134  (
135  "C2",
136  this->coeffDict_,
137  1.92
138  )
139  ),
140  C3_
141  (
142  dimensioned<scalar>::lookupOrAddToDict
143  (
144  "C3",
145  this->coeffDict_,
146  0
147  )
148  ),
149  sigmak_
150  (
151  dimensioned<scalar>::lookupOrAddToDict
152  (
153  "sigmak",
154  this->coeffDict_,
155  1.0
156  )
157  ),
158  sigmaEps_
159  (
160  dimensioned<scalar>::lookupOrAddToDict
161  (
162  "sigmaEps",
163  this->coeffDict_,
164  1.3
165  )
166  ),
167 
168  k_
169  (
170  IOobject
171  (
172  this->groupName("k"),
173  this->runTime_.name(),
174  this->mesh_,
175  IOobject::MUST_READ,
176  IOobject::AUTO_WRITE
177  ),
178  this->mesh_
179  ),
180  epsilon_
181  (
182  IOobject
183  (
184  this->groupName("epsilon"),
185  this->runTime_.name(),
186  this->mesh_,
187  IOobject::MUST_READ,
188  IOobject::AUTO_WRITE
189  ),
190  this->mesh_
191  )
192 {
193  bound(k_, this->kMin_);
194  boundEpsilon();
195 
196  if (type == typeName)
197  {
198  this->printCoeffs(type);
199  }
200 }
201 
202 
203 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
204 
205 template<class BasicMomentumTransportModel>
207 {
209  {
210  Cmu_.readIfPresent(this->coeffDict());
211  C1_.readIfPresent(this->coeffDict());
212  C2_.readIfPresent(this->coeffDict());
213  C3_.readIfPresent(this->coeffDict());
214  sigmak_.readIfPresent(this->coeffDict());
215  sigmaEps_.readIfPresent(this->coeffDict());
216 
217  return true;
218  }
219  else
220  {
221  return false;
222  }
223 }
224 
225 
226 template<class BasicMomentumTransportModel>
228 {
229  if (!this->turbulence_)
230  {
231  return;
232  }
233 
234  // Local references
235  const alphaField& alpha = this->alpha_;
236  const rhoField& rho = this->rho_;
237  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
238  const volVectorField& U = this->U_;
239  volScalarField& nut = this->nut_;
240  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
242  (
243  Foam::fvConstraints::New(this->mesh_)
244  );
245 
247 
249  (
250  fvc::div(fvc::absolute(this->phi(), U))()
251  );
252 
253  tmp<volTensorField> tgradU = fvc::grad(U);
255  (
256  this->GName(),
257  nut()*(dev(twoSymm(tgradU().v())) && tgradU().v())
258  );
259  tgradU.clear();
260 
261  // Update epsilon and G at the wall
262  epsilon_.boundaryFieldRef().updateCoeffs();
263 
264  // Dissipation equation
265  tmp<fvScalarMatrix> epsEqn
266  (
267  fvm::ddt(alpha, rho, epsilon_)
268  + fvm::div(alphaRhoPhi, epsilon_)
269  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
270  ==
271  C1_*alpha()*rho()*G*epsilon_()/k_()
272  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
273  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
274  + epsilonSource()
275  + fvModels.source(alpha, rho, epsilon_)
276  );
277 
278  epsEqn.ref().relax();
279  fvConstraints.constrain(epsEqn.ref());
280  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
281  solve(epsEqn);
282  fvConstraints.constrain(epsilon_);
283  boundEpsilon();
284 
285  // Turbulent kinetic energy equation
287  (
288  fvm::ddt(alpha, rho, k_)
289  + fvm::div(alphaRhoPhi, k_)
290  - fvm::laplacian(alpha*rho*DkEff(), k_)
291  ==
292  alpha()*rho()*G
293  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
294  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
295  + kSource()
296  + fvModels.source(alpha, rho, k_)
297  );
298 
299  kEqn.ref().relax();
300  fvConstraints.constrain(kEqn.ref());
301  solve(kEqn);
303  bound(k_, this->kMin_);
304 
305  correctNut();
306 }
307 
308 
309 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
310 
311 } // End namespace RASModels
312 } // End namespace Foam
313 
314 // ************************************************************************* //
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: kEpsilon.C:74
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:227
tmp< volScalarField > boundEpsilon()
Bound epsilon and return Cmu*sqr(k) for nut.
Definition: kEpsilon.C:41
kEpsilon(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: kEpsilon.C:92
virtual void correctNut()
Correct the eddy-viscosity nut.
Definition: kEpsilon.C:50
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
Definition: kEpsilon.C:59
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilon.C:206
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
const dimensionSet dimVolume
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
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.