LaunderSharmaKE.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-2018 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 "LaunderSharmaKE.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 BasicTurbulenceModel>
41 {
42  return exp(-3.4/sqr(scalar(1) + sqr(k_)/(this->nu()*epsilon_)/50.0));
43 }
44 
45 
46 template<class BasicTurbulenceModel>
48 {
49  return
50  scalar(1)
51  - 0.3*exp(-min(sqr(sqr(k_)/(this->nu()*epsilon_)), scalar(50.0)));
52 }
53 
54 
55 template<class BasicTurbulenceModel>
57 {
58  this->nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
59  this->nut_.correctBoundaryConditions();
60  fv::options::New(this->mesh_).correct(this->nut_);
61 
62  BasicTurbulenceModel::correctNut();
63 }
64 
65 
66 template<class BasicTurbulenceModel>
68 {
69  return tmp<fvScalarMatrix>
70  (
71  new fvScalarMatrix
72  (
73  k_,
74  dimVolume*this->rho_.dimensions()*k_.dimensions()
75  /dimTime
76  )
77  );
78 }
79 
80 
81 template<class BasicTurbulenceModel>
83 {
84  return tmp<fvScalarMatrix>
85  (
86  new fvScalarMatrix
87  (
88  epsilon_,
89  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
90  /dimTime
91  )
92  );
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
97 
98 template<class BasicTurbulenceModel>
100 (
101  const alphaField& alpha,
102  const rhoField& rho,
103  const volVectorField& U,
104  const surfaceScalarField& alphaRhoPhi,
105  const surfaceScalarField& phi,
106  const transportModel& transport,
107  const word& propertiesName,
108  const word& type
109 )
110 :
112  (
113  type,
114  alpha,
115  rho,
116  U,
117  alphaRhoPhi,
118  phi,
119  transport,
120  propertiesName
121  ),
122 
123  Cmu_
124  (
126  (
127  "Cmu",
128  this->coeffDict_,
129  0.09
130  )
131  ),
132  C1_
133  (
135  (
136  "C1",
137  this->coeffDict_,
138  1.44
139  )
140  ),
141  C2_
142  (
144  (
145  "C2",
146  this->coeffDict_,
147  1.92
148  )
149  ),
150  C3_
151  (
153  (
154  "C3",
155  this->coeffDict_,
156  0
157  )
158  ),
159  sigmak_
160  (
162  (
163  "sigmak",
164  this->coeffDict_,
165  1.0
166  )
167  ),
168  sigmaEps_
169  (
171  (
172  "sigmaEps",
173  this->coeffDict_,
174  1.3
175  )
176  ),
177 
178  k_
179  (
180  IOobject
181  (
182  "k",
183  this->runTime_.timeName(),
184  this->mesh_,
187  ),
188  this->mesh_
189  ),
190 
191  epsilon_
192  (
193  IOobject
194  (
195  "epsilon",
196  this->runTime_.timeName(),
197  this->mesh_,
200  ),
201  this->mesh_
202  )
203 {
204  bound(k_, this->kMin_);
205  bound(epsilon_, this->epsilonMin_);
206 
207  if (type == typeName)
208  {
209  this->printCoeffs(type);
210  }
211 }
212 
213 
214 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
215 
216 template<class BasicTurbulenceModel>
218 {
220  {
221  Cmu_.readIfPresent(this->coeffDict());
222  C1_.readIfPresent(this->coeffDict());
223  C2_.readIfPresent(this->coeffDict());
224  C3_.readIfPresent(this->coeffDict());
225  sigmak_.readIfPresent(this->coeffDict());
226  sigmaEps_.readIfPresent(this->coeffDict());
227 
228  return true;
229  }
230  else
231  {
232  return false;
233  }
234 }
235 
236 
237 template<class BasicTurbulenceModel>
239 {
240  if (!this->turbulence_)
241  {
242  return;
243  }
244 
245  // Local references
246  const alphaField& alpha = this->alpha_;
247  const rhoField& rho = this->rho_;
248  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
249  const volVectorField& U = this->U_;
250  volScalarField& nut = this->nut_;
251  fv::options& fvOptions(fv::options::New(this->mesh_));
252 
254 
256 
257  // Calculate parameters and coefficients for Launder-Sharma low-Reynolds
258  // number model
259 
260  volScalarField E(2.0*this->nu()*nut*fvc::magSqrGradGrad(U));
261  volScalarField D(2.0*this->nu()*magSqr(fvc::grad(sqrt(k_))));
262 
263  tmp<volTensorField> tgradU = fvc::grad(U);
264  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
265  tgradU.clear();
266 
267 
268  // Dissipation equation
269  tmp<fvScalarMatrix> epsEqn
270  (
271  fvm::ddt(alpha, rho, epsilon_)
272  + fvm::div(alphaRhoPhi, epsilon_)
273  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
274  ==
275  C1_*alpha*rho*G*epsilon_/k_
276  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha*rho*divU, epsilon_)
277  - fvm::Sp(C2_*f2()*alpha*rho*epsilon_/k_, epsilon_)
278  + alpha*rho*E
279  + epsilonSource()
280  + fvOptions(alpha, rho, epsilon_)
281  );
282 
283  epsEqn.ref().relax();
284  fvOptions.constrain(epsEqn.ref());
285  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
286  solve(epsEqn);
287  fvOptions.correct(epsilon_);
288  bound(epsilon_, this->epsilonMin_);
289 
290 
291  // Turbulent kinetic energy equation
293  (
294  fvm::ddt(alpha, rho, k_)
295  + fvm::div(alphaRhoPhi, k_)
296  - fvm::laplacian(alpha*rho*DkEff(), k_)
297  ==
298  alpha*rho*G - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
299  - fvm::Sp(alpha*rho*(epsilon_ + D)/k_, k_)
300  + kSource()
301  + fvOptions(alpha, rho, k_)
302  );
303 
304  kEqn.ref().relax();
305  fvOptions.constrain(kEqn.ref());
306  solve(kEqn);
307  fvOptions.correct(k_);
308  bound(k_, this->kMin_);
309 
310  correctNut();
311 }
312 
313 
314 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
315 
316 } // End namespace RASModels
317 } // End namespace Foam
318 
319 // ************************************************************************* //
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
surfaceScalarField & phi
virtual bool read()
Re-read model coefficients if they have changed.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
const dimensionedScalar G
Newtonian constant of gravitation.
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
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.
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
virtual tmp< fvScalarMatrix > kSource() const
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
bool read(const char *, int32_t &)
Definition: int32IO.C:85
virtual tmp< fvScalarMatrix > epsilonSource() const
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
A class for handling words, derived from string.
Definition: word.H:59
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
tmp< volScalarField > f2() const
void constrain(fvMatrix< Type > &eqn)
Apply constraints to equation.
dimensioned< scalar > magSqr(const dimensioned< Type > &)
BasicTurbulenceModel::transportModel transportModel
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(psi *p+alphal *rhol0+((alphav *psiv+alphal *psil) - psi) *pSat, rhoMin);# 1 "/home/ubuntu/OpenFOAM-6/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:72
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
fileName::Type type(const fileName &, const bool followLink=true)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:481
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
BasicTurbulenceModel::rhoField rhoField
BasicTurbulenceModel::alphaField alphaField
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
zeroField divU
Definition: alphaSuSp.H:3
tmp< volScalarField > fMu() const
A class for managing temporary objects.
Definition: PtrList.H:53
static options & New(const fvMesh &mesh)
Construct fvOptions and register to datbase if not present.
Definition: fvOptions.C:101
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
Launder and Sharma low-Reynolds k-epsilon turbulence model for incompressible and compressible and co...
volScalarField & nu
scalar nut
Namespace for OpenFOAM.