LaunderSharmaKE.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 "LaunderSharmaKE.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  return exp(-3.4/sqr(scalar(1) + sqr(k_)/(this->nu()*epsilon_)/50.0));
46 }
47 
48 
49 template<class BasicTurbulenceModel>
51 {
52  return
53  scalar(1)
54  - 0.3*exp(-min(sqr(sqr(k_)/(this->nu()*epsilon_)), scalar(50)));
55 }
56 
57 
58 template<class BasicTurbulenceModel>
60 {
61  this->nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
62  this->nut_.correctBoundaryConditions();
63  fv::options::New(this->mesh_).correct(this->nut_);
64 
65  BasicTurbulenceModel::correctNut();
66 }
67 
68 
69 template<class BasicTurbulenceModel>
71 {
73  (
74  k_,
75  dimVolume*this->rho_.dimensions()*k_.dimensions()/dimTime
76  );
77 }
78 
79 
80 template<class BasicTurbulenceModel>
82 {
84  (
85  epsilon_,
86  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()/dimTime
87  );
88 }
89 
90 
91 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
92 
93 template<class BasicTurbulenceModel>
95 (
96  const alphaField& alpha,
97  const rhoField& rho,
98  const volVectorField& U,
99  const surfaceScalarField& alphaRhoPhi,
100  const surfaceScalarField& phi,
101  const transportModel& transport,
102  const word& propertiesName,
103  const word& type
104 )
105 :
106  eddyViscosity<RASModel<BasicTurbulenceModel>>
107  (
108  type,
109  alpha,
110  rho,
111  U,
112  alphaRhoPhi,
113  phi,
114  transport,
115  propertiesName
116  ),
117 
118  Cmu_
119  (
120  dimensioned<scalar>::getOrAddToDict
121  (
122  "Cmu",
123  this->coeffDict_,
124  0.09
125  )
126  ),
127  C1_
128  (
129  dimensioned<scalar>::getOrAddToDict
130  (
131  "C1",
132  this->coeffDict_,
133  1.44
134  )
135  ),
136  C2_
137  (
138  dimensioned<scalar>::getOrAddToDict
139  (
140  "C2",
141  this->coeffDict_,
142  1.92
143  )
144  ),
145  C3_
146  (
147  dimensioned<scalar>::getOrAddToDict
148  (
149  "C3",
150  this->coeffDict_,
151  0
152  )
153  ),
154  sigmak_
155  (
156  dimensioned<scalar>::getOrAddToDict
157  (
158  "sigmak",
159  this->coeffDict_,
160  1.0
161  )
162  ),
163  sigmaEps_
164  (
165  dimensioned<scalar>::getOrAddToDict
166  (
167  "sigmaEps",
168  this->coeffDict_,
169  1.3
170  )
171  ),
172 
173  k_
174  (
175  IOobject
176  (
177  "k",
178  this->runTime_.timeName(),
179  this->mesh_,
180  IOobject::MUST_READ,
181  IOobject::AUTO_WRITE
182  ),
183  this->mesh_
184  ),
185 
186  epsilon_
187  (
188  IOobject
189  (
190  "epsilon",
191  this->runTime_.timeName(),
192  this->mesh_,
193  IOobject::MUST_READ,
194  IOobject::AUTO_WRITE
195  ),
196  this->mesh_
197  )
198 {
199  bound(k_, this->kMin_);
200  bound(epsilon_, this->epsilonMin_);
201 
202  if (type == typeName)
203  {
204  this->printCoeffs(type);
205  }
206 }
207 
208 
209 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
210 
211 template<class BasicTurbulenceModel>
213 {
215  {
216  Cmu_.readIfPresent(this->coeffDict());
217  C1_.readIfPresent(this->coeffDict());
218  C2_.readIfPresent(this->coeffDict());
219  C3_.readIfPresent(this->coeffDict());
220  sigmak_.readIfPresent(this->coeffDict());
221  sigmaEps_.readIfPresent(this->coeffDict());
222 
223  return true;
224  }
225 
226  return false;
227 }
228 
229 
230 template<class BasicTurbulenceModel>
232 {
233  if (!this->turbulence_)
234  {
235  return;
236  }
237 
238  // Local references
239  const alphaField& alpha = this->alpha_;
240  const rhoField& rho = this->rho_;
241  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
242  const volVectorField& U = this->U_;
243  volScalarField& nut = this->nut_;
244  fv::options& fvOptions(fv::options::New(this->mesh_));
245 
246  eddyViscosity<RASModel<BasicTurbulenceModel>>::correct();
247 
249 
250  // Calculate parameters and coefficients for Launder-Sharma low-Reynolds
251  // number model
252 
253  volScalarField E(2.0*this->nu()*nut*fvc::magSqrGradGrad(U));
254  volScalarField D(2.0*this->nu()*magSqr(fvc::grad(sqrt(k_))));
255 
256  tmp<volTensorField> tgradU = fvc::grad(U);
257  volScalarField G(this->GName(), nut*(tgradU() && devTwoSymm(tgradU())));
258  tgradU.clear();
259 
260 
261  // Dissipation equation
262  tmp<fvScalarMatrix> epsEqn
263  (
264  fvm::ddt(alpha, rho, epsilon_)
265  + fvm::div(alphaRhoPhi, epsilon_)
266  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
267  ==
268  C1_*alpha*rho*G*epsilon_/k_
269  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha*rho*divU, epsilon_)
270  - fvm::Sp(C2_*f2()*alpha*rho*epsilon_/k_, epsilon_)
271  + alpha*rho*E
272  + epsilonSource()
273  + fvOptions(alpha, rho, epsilon_)
274  );
275 
276  epsEqn.ref().relax();
277  fvOptions.constrain(epsEqn.ref());
278  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
279  solve(epsEqn);
280  fvOptions.correct(epsilon_);
281  bound(epsilon_, this->epsilonMin_);
282 
283 
284  // Turbulent kinetic energy equation
285  tmp<fvScalarMatrix> kEqn
286  (
287  fvm::ddt(alpha, rho, k_)
288  + fvm::div(alphaRhoPhi, k_)
289  - fvm::laplacian(alpha*rho*DkEff(), k_)
290  ==
291  alpha*rho*G - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
292  - fvm::Sp(alpha*rho*(epsilon_ + D)/k_, k_)
293  + kSource()
294  + fvOptions(alpha, rho, k_)
295  );
296 
297  kEqn.ref().relax();
298  fvOptions.constrain(kEqn.ref());
299  solve(kEqn);
300  fvOptions.correct(k_);
301  bound(k_, this->kMin_);
302 
303  correctNut();
304 }
305 
306 
307 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
308 
309 } // End namespace RASModels
310 } // End namespace Foam
311 
312 // ************************************************************************* //
void correct(GeometricField< Type, PatchField, GeoMesh > &field)
Apply correction to field.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
dimensionedScalar kMin_
Lower limit of k.
Definition: RASModel.H:109
virtual bool read()
Re-read model coefficients if they have changed.
const dimensionedScalar G
Newtonian constant of gravitation.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< volScalarField > magSqrGradGrad(const GeometricField< Type, fvPatchField, volMesh > &vf)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:51
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:76
virtual tmp< fvScalarMatrix > kSource() const
fv::options & fvOptions
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:77
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
virtual tmp< fvScalarMatrix > epsilonSource() const
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from Foam::string.
Definition: word.H:63
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
tmp< volScalarField > f2() const
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:27
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
SolverPerformance< Type > solve(faMatrix< Type > &, const dictionary &solverControls)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
Bound the given scalar field if it has gone unbounded.
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
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:41
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
const dimensionedScalar & D
zeroField divU
Definition: alphaSuSp.H:3
tmp< volScalarField > fMu() const
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
Launder and Sharma low-Reynolds k-epsilon turbulence model for incompressible and compressible and co...
volScalarField & nu
dimensionedScalar epsilonMin_
Lower limit of epsilon.
Definition: RASModel.H:114
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
scalar nut
Namespace for OpenFOAM.
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491