kEpsilon.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 "kEpsilon.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()
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 :
98  eddyViscosity<RASModel<BasicTurbulenceModel>>
99  (
100  type,
101  alpha,
102  rho,
103  U,
104  alphaRhoPhi,
105  phi,
106  transport,
107  propertiesName
108  ),
109 
110  Cmu_
111  (
112  dimensioned<scalar>::getOrAddToDict
113  (
114  "Cmu",
115  this->coeffDict_,
116  0.09
117  )
118  ),
119  C1_
120  (
121  dimensioned<scalar>::getOrAddToDict
122  (
123  "C1",
124  this->coeffDict_,
125  1.44
126  )
127  ),
128  C2_
129  (
130  dimensioned<scalar>::getOrAddToDict
131  (
132  "C2",
133  this->coeffDict_,
134  1.92
135  )
136  ),
137  C3_
138  (
139  dimensioned<scalar>::getOrAddToDict
140  (
141  "C3",
142  this->coeffDict_,
143  0
144  )
145  ),
146  sigmak_
147  (
148  dimensioned<scalar>::getOrAddToDict
149  (
150  "sigmak",
151  this->coeffDict_,
152  1.0
153  )
154  ),
155  sigmaEps_
156  (
157  dimensioned<scalar>::getOrAddToDict
158  (
159  "sigmaEps",
160  this->coeffDict_,
161  1.3
162  )
163  ),
164 
165  k_
166  (
167  IOobject
168  (
169  IOobject::groupName("k", alphaRhoPhi.group()),
170  this->runTime_.timeName(),
171  this->mesh_,
172  IOobject::MUST_READ,
173  IOobject::AUTO_WRITE
174  ),
175  this->mesh_
176  ),
177  epsilon_
178  (
179  IOobject
180  (
181  IOobject::groupName("epsilon", alphaRhoPhi.group()),
182  this->runTime_.timeName(),
183  this->mesh_,
184  IOobject::MUST_READ,
185  IOobject::AUTO_WRITE
186  ),
187  this->mesh_
188  )
189 {
190  bound(k_, this->kMin_);
191  bound(epsilon_, this->epsilonMin_);
192 
193  if (type == typeName)
194  {
195  this->printCoeffs(type);
196  }
197 }
198 
199 
200 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
201 
202 template<class BasicTurbulenceModel>
204 {
206  {
207  Cmu_.readIfPresent(this->coeffDict());
208  C1_.readIfPresent(this->coeffDict());
209  C2_.readIfPresent(this->coeffDict());
210  C3_.readIfPresent(this->coeffDict());
211  sigmak_.readIfPresent(this->coeffDict());
212  sigmaEps_.readIfPresent(this->coeffDict());
213 
214  return true;
215  }
216 
217  return false;
218 }
219 
220 
221 template<class BasicTurbulenceModel>
223 {
224  if (!this->turbulence_)
225  {
226  return;
227  }
228 
229  // Local references
230  const alphaField& alpha = this->alpha_;
231  const rhoField& rho = this->rho_;
232  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
233  const volVectorField& U = this->U_;
234  const volScalarField& nut = this->nut_;
235 
236  fv::options& fvOptions(fv::options::New(this->mesh_));
237 
238  eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
239 
240  const volScalarField::Internal divU
241  (
242  fvc::div(fvc::absolute(this->phi(), U))().v()
243  );
244 
245  tmp<volTensorField> tgradU = fvc::grad(U);
246  const volScalarField::Internal GbyNu
247  (
248  IOobject::scopedName(this->type(), "GbyNu"),
249  tgradU().v() && devTwoSymm(tgradU().v())
250  );
251  const volScalarField::Internal G(this->GName(), nut()*GbyNu);
252  tgradU.clear();
253 
254  // Update epsilon and G at the wall
255  epsilon_.boundaryFieldRef().updateCoeffs();
256  // Push any changed cell values to coupled neighbours
257  epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
258 
259  // Dissipation equation
260  tmp<fvScalarMatrix> epsEqn
261  (
262  fvm::ddt(alpha, rho, epsilon_)
263  + fvm::div(alphaRhoPhi, epsilon_)
264  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
265  ==
266  C1_*alpha()*rho()*GbyNu*Cmu_*k_()
267  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
268  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
269  + epsilonSource()
270  + fvOptions(alpha, rho, epsilon_)
271  );
272 
273  epsEqn.ref().relax();
274  fvOptions.constrain(epsEqn.ref());
275  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
276  solve(epsEqn);
277  fvOptions.correct(epsilon_);
278  bound(epsilon_, this->epsilonMin_);
279 
280  // Turbulent kinetic energy equation
281  tmp<fvScalarMatrix> kEqn
282  (
283  fvm::ddt(alpha, rho, k_)
284  + fvm::div(alphaRhoPhi, k_)
285  - fvm::laplacian(alpha*rho*DkEff(), k_)
286  ==
287  alpha()*rho()*G
288  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
289  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
290  + kSource()
291  + fvOptions(alpha, rho, k_)
292  );
293 
294  kEqn.ref().relax();
295  fvOptions.constrain(kEqn.ref());
296  solve(kEqn);
297  fvOptions.correct(k_);
298  bound(k_, this->kMin_);
299 
300  correctNut();
301 }
302 
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 } // End namespace RASModels
307 } // End namespace Foam
308 
309 // ************************************************************************* //
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:62
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
type
Types of root.
Definition: Roots.H:52
const dimensionedScalar G
Newtonian constant of gravitation.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, faPatchField, areaMesh > > div(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facDiv.C:43
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:47
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
Generic dimensioned Type class.
tmp< GeometricField< Type, faPatchField, areaMesh > > ddt(const dimensioned< Type > dt, const faMesh &mesh)
Definition: facDdt.C:40
zeroField SuSp(const Foam::zero, const GeometricField< Type, faPatchField, areaMesh > &)
A no-op source.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
fv::options & fvOptions
CEqn solve()
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:81
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
Standard k-epsilon turbulence model for incompressible and compressible flows including rapid distort...
Definition: kEpsilon.H:85
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< GeometricField< Type, faPatchField, areaMesh > > laplacian(const GeometricField< Type, faPatchField, areaMesh > &vf, const word &name)
Definition: facLaplacian.C:40
Bound the given scalar field if it has gone unbounded.
void read(Istream &, label &val, const dictionary &)
In-place read with dictionary lookup.
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/v2312/OpenFOAM-v2312/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
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
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
virtual void correctNut()
Definition: kEpsilon.C:36
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:172
scalar nut
Namespace for OpenFOAM.
tmp< GeometricField< typename outerProduct< vector, Type >::type, faPatchField, areaMesh >> grad(const GeometricField< Type, faePatchField, edgeMesh > &ssf)
Definition: facGrad.C:51
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491
zeroField Sp
Definition: alphaSuSp.H:2