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 | Copyright (C) 2011-2015 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 "bound.H"
28 
29 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace RASModels
34 {
35 
36 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
37 
38 template<class BasicTurbulenceModel>
40 {
41  return exp(-3.4/sqr(scalar(1) + sqr(k_)/(this->nu()*epsilon_)/50.0));
42 }
43 
44 
45 template<class BasicTurbulenceModel>
47 {
48  return
49  scalar(1)
50  - 0.3*exp(-min(sqr(sqr(k_)/(this->nu()*epsilon_)), scalar(50.0)));
51 }
52 
53 
54 template<class BasicTurbulenceModel>
56 {
57  this->nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
58  this->nut_.correctBoundaryConditions();
59 
60  BasicTurbulenceModel::correctNut();
61 }
62 
63 
64 template<class BasicTurbulenceModel>
66 {
67  return tmp<fvScalarMatrix>
68  (
69  new fvScalarMatrix
70  (
71  k_,
72  dimVolume*this->rho_.dimensions()*k_.dimensions()
73  /dimTime
74  )
75  );
76 }
77 
78 
79 template<class BasicTurbulenceModel>
81 {
82  return tmp<fvScalarMatrix>
83  (
84  new fvScalarMatrix
85  (
86  epsilon_,
87  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
88  /dimTime
89  )
90  );
91 }
92 
93 
94 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
95 
96 template<class BasicTurbulenceModel>
98 (
99  const alphaField& alpha,
100  const rhoField& rho,
101  const volVectorField& U,
102  const surfaceScalarField& alphaRhoPhi,
103  const surfaceScalarField& phi,
104  const transportModel& transport,
105  const word& propertiesName,
106  const word& type
107 )
108 :
110  (
111  type,
112  alpha,
113  rho,
114  U,
115  alphaRhoPhi,
116  phi,
117  transport,
118  propertiesName
119  ),
120 
121  Cmu_
122  (
124  (
125  "Cmu",
126  this->coeffDict_,
127  0.09
128  )
129  ),
130  C1_
131  (
133  (
134  "C1",
135  this->coeffDict_,
136  1.44
137  )
138  ),
139  C2_
140  (
142  (
143  "C2",
144  this->coeffDict_,
145  1.92
146  )
147  ),
148  C3_
149  (
151  (
152  "C3",
153  this->coeffDict_,
154  -0.33
155  )
156  ),
157  sigmak_
158  (
160  (
161  "sigmak",
162  this->coeffDict_,
163  1.0
164  )
165  ),
166  sigmaEps_
167  (
169  (
170  "sigmaEps",
171  this->coeffDict_,
172  1.3
173  )
174  ),
175 
176  k_
177  (
178  IOobject
179  (
180  "k",
181  this->runTime_.timeName(),
182  this->mesh_,
185  ),
186  this->mesh_
187  ),
188 
189  epsilon_
190  (
191  IOobject
192  (
193  "epsilon",
194  this->runTime_.timeName(),
195  this->mesh_,
198  ),
199  this->mesh_
200  )
201 {
202  bound(k_, this->kMin_);
203  bound(epsilon_, this->epsilonMin_);
204 
205  if (type == typeName)
206  {
207  correctNut();
208  this->printCoeffs(type);
209  }
210 }
211 
212 
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
214 
215 template<class BasicTurbulenceModel>
217 {
219  {
220  Cmu_.readIfPresent(this->coeffDict());
221  C1_.readIfPresent(this->coeffDict());
222  C2_.readIfPresent(this->coeffDict());
223  C3_.readIfPresent(this->coeffDict());
224  sigmak_.readIfPresent(this->coeffDict());
225  sigmaEps_.readIfPresent(this->coeffDict());
226 
227  return true;
228  }
229  else
230  {
231  return false;
232  }
233 }
234 
235 
236 template<class BasicTurbulenceModel>
238 {
239  if (!this->turbulence_)
240  {
241  return;
242  }
243 
244  // Local references
245  const alphaField& alpha = this->alpha_;
246  const rhoField& rho = this->rho_;
247  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
248  const volVectorField& U = this->U_;
249  volScalarField& nut = this->nut_;
250 
252 
254 
255  // Calculate parameters and coefficients for Launder-Sharma low-Reynolds
256  // number model
257 
258  volScalarField E(2.0*this->nu()*nut*fvc::magSqrGradGrad(U));
259  volScalarField D(2.0*this->nu()*magSqr(fvc::grad(sqrt(k_))));
260 
261  tmp<volTensorField> tgradU = fvc::grad(U);
262  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
263  tgradU.clear();
264 
265 
266  // Dissipation equation
267  tmp<fvScalarMatrix> epsEqn
268  (
269  fvm::ddt(alpha, rho, epsilon_)
270  + fvm::div(alphaRhoPhi, epsilon_)
271  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
272  ==
273  C1_*alpha*rho*G*epsilon_/k_
274  - fvm::SuSp(((2.0/3.0)*C1_ + C3_)*alpha*rho*divU, epsilon_)
275  - fvm::Sp(C2_*f2()*alpha*rho*epsilon_/k_, epsilon_)
276  + alpha*rho*E
277  + epsilonSource()
278  );
279 
280  epsEqn().relax();
281  solve(epsEqn);
282  bound(epsilon_, this->epsilonMin_);
283 
284 
285  // Turbulent kinetic energy equation
287  (
288  fvm::ddt(alpha, rho, k_)
289  + fvm::div(alphaRhoPhi, k_)
290  - fvm::laplacian(alpha*rho*DkEff(), k_)
291  ==
292  alpha*rho*G - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
293  - fvm::Sp(alpha*rho*(epsilon_ + D)/k_, k_)
294  + kSource()
295  );
296 
297  kEqn().relax();
298  solve(kEqn);
299  bound(k_, this->kMin_);
300 
301  correctNut();
302 }
303 
304 
305 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
306 
307 } // End namespace RASModels
308 } // End namespace Foam
309 
310 // ************************************************************************* //
dimensionedScalar sqrt(const dimensionedScalar &ds)
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
BasicTurbulenceModel::transportModel transportModel
BasicTurbulenceModel::alphaField alphaField
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
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.
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:174
dimensioned< scalar > magSqr(const dimensioned< Type > &)
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
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/dm2/henry/OpenFOAM/OpenFOAM-3.0.x/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:74
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
tmp< volScalarField > fMu() const
A class for handling words, derived from string.
Definition: word.H:59
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
tmp< volScalarField > magSqrGradGrad(const GeometricField< Type, fvPatchField, volMesh > &vf)
dimensionedScalar exp(const dimensionedScalar &ds)
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
virtual tmp< fvScalarMatrix > kSource() const
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Generic dimensioned Type class.
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
scalar nut
tmp< volScalarField > f2() const
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:589
surfaceScalarField & phi
virtual bool read()
Re-read model coefficients if they have changed.
solverPerformance solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
virtual tmp< fvScalarMatrix > epsilonSource() const
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
bool read(const char *, int32_t &)
Definition: int32IO.C:87
volScalarField & nu
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
BasicTurbulenceModel::rhoField rhoField
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
A class for managing temporary objects.
Definition: PtrList.H:118
const dimensionedScalar G
Newtonian constant of gravitation.
U
Definition: pEqn.H:82
Launder and Sharma low-Reynolds k-epsilon turbulence model for incompressible and compressible and co...