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-2020 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 "fvOptions.H"
28 #include "bound.H"
29 
30 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace RASModels
35 {
36 
37 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
38 
39 template<class BasicMomentumTransportModel>
41 {
42  this->nut_ = Cmu_*sqr(k_)/epsilon_;
43  this->nut_.correctBoundaryConditions();
44  fv::options::New(this->mesh_).correct(this->nut_);
45 }
46 
47 
48 template<class BasicMomentumTransportModel>
50 {
51  return tmp<fvScalarMatrix>
52  (
53  new fvScalarMatrix
54  (
55  k_,
56  dimVolume*this->rho_.dimensions()*k_.dimensions()
57  /dimTime
58  )
59  );
60 }
61 
62 
63 template<class BasicMomentumTransportModel>
65 {
66  return tmp<fvScalarMatrix>
67  (
68  new fvScalarMatrix
69  (
70  epsilon_,
71  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
72  /dimTime
73  )
74  );
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
79 
80 template<class BasicMomentumTransportModel>
82 (
83  const alphaField& alpha,
84  const rhoField& rho,
85  const volVectorField& U,
86  const surfaceScalarField& alphaRhoPhi,
87  const surfaceScalarField& phi,
88  const transportModel& transport,
89  const word& type
90 )
91 :
93  (
94  type,
95  alpha,
96  rho,
97  U,
98  alphaRhoPhi,
99  phi,
100  transport
101  ),
102 
103  Cmu_
104  (
106  (
107  "Cmu",
108  this->coeffDict_,
109  0.09
110  )
111  ),
112  C1_
113  (
115  (
116  "C1",
117  this->coeffDict_,
118  1.44
119  )
120  ),
121  C2_
122  (
124  (
125  "C2",
126  this->coeffDict_,
127  1.92
128  )
129  ),
130  C3_
131  (
133  (
134  "C3",
135  this->coeffDict_,
136  0
137  )
138  ),
139  sigmak_
140  (
142  (
143  "sigmak",
144  this->coeffDict_,
145  1.0
146  )
147  ),
148  sigmaEps_
149  (
151  (
152  "sigmaEps",
153  this->coeffDict_,
154  1.3
155  )
156  ),
157 
158  k_
159  (
160  IOobject
161  (
162  IOobject::groupName("k", alphaRhoPhi.group()),
163  this->runTime_.timeName(),
164  this->mesh_,
167  ),
168  this->mesh_
169  ),
170  epsilon_
171  (
172  IOobject
173  (
174  IOobject::groupName("epsilon", alphaRhoPhi.group()),
175  this->runTime_.timeName(),
176  this->mesh_,
179  ),
180  this->mesh_
181  )
182 {
183  bound(k_, this->kMin_);
184  bound(epsilon_, this->epsilonMin_);
185 
186  if (type == typeName)
187  {
188  this->printCoeffs(type);
189  }
190 }
191 
192 
193 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
194 
195 template<class BasicMomentumTransportModel>
197 {
199  {
200  Cmu_.readIfPresent(this->coeffDict());
201  C1_.readIfPresent(this->coeffDict());
202  C2_.readIfPresent(this->coeffDict());
203  C3_.readIfPresent(this->coeffDict());
204  sigmak_.readIfPresent(this->coeffDict());
205  sigmaEps_.readIfPresent(this->coeffDict());
206 
207  return true;
208  }
209  else
210  {
211  return false;
212  }
213 }
214 
215 
216 template<class BasicMomentumTransportModel>
218 {
219  if (!this->turbulence_)
220  {
221  return;
222  }
223 
224  // Local references
225  const alphaField& alpha = this->alpha_;
226  const rhoField& rho = this->rho_;
227  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
228  const volVectorField& U = this->U_;
229  volScalarField& nut = this->nut_;
230  fv::options& fvOptions(fv::options::New(this->mesh_));
231 
233 
235  (
236  fvc::div(fvc::absolute(this->phi(), U))()
237  );
238 
239  tmp<volTensorField> tgradU = fvc::grad(U);
241  (
242  this->GName(),
243  nut()*(dev(twoSymm(tgradU().v())) && tgradU().v())
244  );
245  tgradU.clear();
246 
247  // Update epsilon and G at the wall
248  epsilon_.boundaryFieldRef().updateCoeffs();
249 
250  // Dissipation equation
251  tmp<fvScalarMatrix> epsEqn
252  (
253  fvm::ddt(alpha, rho, epsilon_)
254  + fvm::div(alphaRhoPhi, epsilon_)
255  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
256  ==
257  C1_*alpha()*rho()*G*epsilon_()/k_()
258  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
259  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
260  + epsilonSource()
261  + fvOptions(alpha, rho, epsilon_)
262  );
263 
264  epsEqn.ref().relax();
265  fvOptions.constrain(epsEqn.ref());
266  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
267  solve(epsEqn);
268  fvOptions.correct(epsilon_);
269  bound(epsilon_, this->epsilonMin_);
270 
271  // Turbulent kinetic energy equation
273  (
274  fvm::ddt(alpha, rho, k_)
275  + fvm::div(alphaRhoPhi, k_)
276  - fvm::laplacian(alpha*rho*DkEff(), k_)
277  ==
278  alpha()*rho()*G
279  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
280  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
281  + kSource()
282  + fvOptions(alpha, rho, k_)
283  );
284 
285  kEqn.ref().relax();
286  fvOptions.constrain(kEqn.ref());
287  solve(kEqn);
288  fvOptions.correct(k_);
289  bound(k_, this->kMin_);
290 
291  correctNut();
292 }
293 
294 
295 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
296 
297 } // End namespace RASModels
298 } // End namespace Foam
299 
300 // ************************************************************************* //
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:64
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:230
fv::options & fvOptions
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:174
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
virtual tmp< fvScalarMatrix > kSource() const
Definition: kEpsilon.C:49
Generic dimensioned Type class.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
Finite-volume options.
Definition: fvOptions.H:52
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
kEpsilon(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
Definition: kEpsilon.C:82
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
bool read(const char *, int32_t &)
Definition: int32IO.C:85
A class representing the concept of a GeometricField of 1 used to avoid unnecessary manipulations for...
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Bound the given scalar field if it has gone unbounded.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
virtual bool read()
Re-read model coefficients if they have changed.
Definition: kEpsilon.C:196
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:188
U
Definition: pEqn.H:72
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
virtual void correctNut()
Definition: kEpsilon.C:40
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
zeroField divU
Definition: alphaSuSp.H:3
A class for managing temporary objects.
Definition: PtrList.H:53
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: kEpsilon.C:217
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionedScalar & G
Newtonian constant of gravitation.
scalar nut
Namespace for OpenFOAM.