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 {
57  (
58  k_,
59  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
60  );
61 }
62 
63 
64 template<class BasicTurbulenceModel>
66 {
68  (
69  epsilon_,
70  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()/dimTime
71  );
72 }
73 
74 
75 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
76 
77 template<class BasicTurbulenceModel>
79 (
80  const alphaField& alpha,
81  const rhoField& rho,
82  const volVectorField& U,
83  const surfaceScalarField& alphaRhoPhi,
84  const surfaceScalarField& phi,
85  const transportModel& transport,
86  const word& propertiesName,
87  const word& type
88 )
89 :
90  eddyViscosity<RASModel<BasicTurbulenceModel>>
91  (
92  type,
93  alpha,
94  rho,
95  U,
96  alphaRhoPhi,
97  phi,
98  transport,
99  propertiesName
100  ),
101 
102  Cmu_
103  (
104  dimensioned<scalar>::getOrAddToDict
105  (
106  "Cmu",
107  this->coeffDict_,
108  0.09
109  )
110  ),
111  C1_
112  (
113  dimensioned<scalar>::getOrAddToDict
114  (
115  "C1",
116  this->coeffDict_,
117  1.44
118  )
119  ),
120  C2_
121  (
122  dimensioned<scalar>::getOrAddToDict
123  (
124  "C2",
125  this->coeffDict_,
126  1.92
127  )
128  ),
129  C3_
130  (
131  dimensioned<scalar>::getOrAddToDict
132  (
133  "C3",
134  this->coeffDict_,
135  0
136  )
137  ),
138  sigmak_
139  (
140  dimensioned<scalar>::getOrAddToDict
141  (
142  "sigmak",
143  this->coeffDict_,
144  1.0
145  )
146  ),
147  sigmaEps_
148  (
149  dimensioned<scalar>::getOrAddToDict
150  (
151  "sigmaEps",
152  this->coeffDict_,
153  1.3
154  )
155  ),
156 
157  k_
158  (
159  IOobject
160  (
161  IOobject::groupName("k", alphaRhoPhi.group()),
162  this->runTime_.timeName(),
163  this->mesh_,
164  IOobject::MUST_READ,
165  IOobject::AUTO_WRITE
166  ),
167  this->mesh_
168  ),
169  epsilon_
170  (
171  IOobject
172  (
173  IOobject::groupName("epsilon", alphaRhoPhi.group()),
174  this->runTime_.timeName(),
175  this->mesh_,
176  IOobject::MUST_READ,
177  IOobject::AUTO_WRITE
178  ),
179  this->mesh_
180  )
181 {
182  bound(k_, this->kMin_);
183  bound(epsilon_, this->epsilonMin_);
184 
185  if (type == typeName)
186  {
187  this->printCoeffs(type);
188  }
189 }
190 
191 
192 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
193 
194 template<class BasicTurbulenceModel>
196 {
198  {
199  Cmu_.readIfPresent(this->coeffDict());
200  C1_.readIfPresent(this->coeffDict());
201  C2_.readIfPresent(this->coeffDict());
202  C3_.readIfPresent(this->coeffDict());
203  sigmak_.readIfPresent(this->coeffDict());
204  sigmaEps_.readIfPresent(this->coeffDict());
205 
206  return true;
207  }
208 
209  return false;
210 }
211 
212 
213 template<class BasicTurbulenceModel>
215 {
216  if (!this->turbulence_)
217  {
218  return;
219  }
220 
221  // Local references
222  const alphaField& alpha = this->alpha_;
223  const rhoField& rho = this->rho_;
224  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
225  const volVectorField& U = this->U_;
226  const volScalarField& nut = this->nut_;
227 
228  fv::options& fvOptions(fv::options::New(this->mesh_));
229 
230  eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
231 
232  const volScalarField::Internal divU
233  (
234  fvc::div(fvc::absolute(this->phi(), U))().v()
235  );
236 
237  tmp<volTensorField> tgradU = fvc::grad(U);
238  const volScalarField::Internal GbyNu
239  (
240  IOobject::scopedName(this->type(), "GbyNu"),
241  tgradU().v() && devTwoSymm(tgradU().v())
242  );
243  const volScalarField::Internal G(this->GName(), nut()*GbyNu);
244  tgradU.clear();
245 
246  // Update epsilon and G at the wall
247  epsilon_.boundaryFieldRef().updateCoeffs();
248  // Push any changed cell values to coupled neighbours
249  epsilon_.boundaryFieldRef().template evaluateCoupled<coupledFvPatch>();
250 
251  // Dissipation equation
252  tmp<fvScalarMatrix> epsEqn
253  (
254  fvm::ddt(alpha, rho, epsilon_)
255  + fvm::div(alphaRhoPhi, epsilon_)
256  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
257  ==
258  C1_*alpha()*rho()*GbyNu*Cmu_*k_()
259  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha()*rho()*divU, epsilon_)
260  - fvm::Sp(C2_*alpha()*rho()*epsilon_()/k_(), epsilon_)
261  + epsilonSource()
262  + fvOptions(alpha, rho, epsilon_)
263  );
264 
265  epsEqn.ref().relax();
266  fvOptions.constrain(epsEqn.ref());
267  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
268  solve(epsEqn);
269  fvOptions.correct(epsilon_);
270  bound(epsilon_, this->epsilonMin_);
271 
272  // Turbulent kinetic energy equation
273  tmp<fvScalarMatrix> kEqn
274  (
275  fvm::ddt(alpha, rho, k_)
276  + fvm::div(alphaRhoPhi, k_)
277  - fvm::laplacian(alpha*rho*DkEff(), k_)
278  ==
279  alpha()*rho()*G
280  - fvm::SuSp((2.0/3.0)*alpha()*rho()*divU, k_)
281  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
282  + kSource()
283  + fvOptions(alpha, rho, k_)
284  );
285 
286  kEqn.ref().relax();
287  fvOptions.constrain(kEqn.ref());
288  solve(kEqn);
289  fvOptions.correct(k_);
290  bound(k_, this->kMin_);
291 
292  correctNut();
293 }
294 
295 
296 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
297 
298 } // End namespace RASModels
299 } // End namespace Foam
300 
301 // ************************************************************************* //
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: kEpsilon.C:58
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
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:76
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:72
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
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
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/v2406/OpenFOAM-v2406/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:180
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