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  this->printCoeffs(type);
208  }
209 }
210 
211 
212 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
213 
214 template<class BasicTurbulenceModel>
216 {
218  {
219  Cmu_.readIfPresent(this->coeffDict());
220  C1_.readIfPresent(this->coeffDict());
221  C2_.readIfPresent(this->coeffDict());
222  C3_.readIfPresent(this->coeffDict());
223  sigmak_.readIfPresent(this->coeffDict());
224  sigmaEps_.readIfPresent(this->coeffDict());
225 
226  return true;
227  }
228  else
229  {
230  return false;
231  }
232 }
233 
234 
235 template<class BasicTurbulenceModel>
237 {
238  if (!this->turbulence_)
239  {
240  return;
241  }
242 
243  // Local references
244  const alphaField& alpha = this->alpha_;
245  const rhoField& rho = this->rho_;
246  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
247  const volVectorField& U = this->U_;
248  volScalarField& nut = this->nut_;
249 
251 
253 
254  // Calculate parameters and coefficients for Launder-Sharma low-Reynolds
255  // number model
256 
257  volScalarField E(2.0*this->nu()*nut*fvc::magSqrGradGrad(U));
258  volScalarField D(2.0*this->nu()*magSqr(fvc::grad(sqrt(k_))));
259 
260  tmp<volTensorField> tgradU = fvc::grad(U);
261  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
262  tgradU.clear();
263 
264 
265  // Dissipation equation
266  tmp<fvScalarMatrix> epsEqn
267  (
268  fvm::ddt(alpha, rho, epsilon_)
269  + fvm::div(alphaRhoPhi, epsilon_)
270  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
271  ==
272  C1_*alpha*rho*G*epsilon_/k_
273  - fvm::SuSp(((2.0/3.0)*C1_ + C3_)*alpha*rho*divU, epsilon_)
274  - fvm::Sp(C2_*f2()*alpha*rho*epsilon_/k_, epsilon_)
275  + alpha*rho*E
276  + epsilonSource()
277  );
278 
279  epsEqn().relax();
280  solve(epsEqn);
281  bound(epsilon_, this->epsilonMin_);
282 
283 
284  // Turbulent kinetic energy equation
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  );
295 
296  kEqn().relax();
297  solve(kEqn);
298  bound(k_, this->kMin_);
299 
300  correctNut();
301 }
302 
303 
304 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
305 
306 } // End namespace RASModels
307 } // End namespace Foam
308 
309 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Foam::IOobject::AUTO_WRITE
@ AUTO_WRITE
Definition: IOobject.H:117
Foam::constant::universal::G
const dimensionedScalar G
Newtonian constant of gravitation.
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
LaunderSharmaKE.H
Foam::RASModel
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
phi
surfaceScalarField & phi
Definition: setRegionFluidFields.H:8
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:216
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::bound
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Foam::RASModels::LaunderSharmaKE::correct
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: LaunderSharmaKE.C:236
Foam::fvc::div
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
Foam::IOobject::MUST_READ
@ MUST_READ
Definition: IOobject.H:108
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::RASModels::LaunderSharmaKE::fMu
tmp< volScalarField > fMu() const
Definition: LaunderSharmaKE.C:39
U
U
Definition: pEqn.H:46
Foam::fvm::Sp
tmp< fvMatrix< Type > > Sp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
nu
volScalarField & nu
Definition: readMechanicalProperties.H:179
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::RASModels::LaunderSharmaKE::epsilonSource
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: LaunderSharmaKE.C:80
nut
nut
Definition: createTDFields.H:71
Foam::RASModels::LaunderSharmaKE::alphaField
BasicTurbulenceModel::alphaField alphaField
Definition: LaunderSharmaKE.H:127
Foam::fvm::laplacian
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
Foam::RASModels::LaunderSharmaKE::read
virtual bool read()
Re-read model coefficients if they have changed.
Definition: LaunderSharmaKE.C:215
bound.H
Bound the given scalar field if it has gone unbounded.
correct
fvOptions correct(rho)
Foam::fvm::SuSp
tmp< fvMatrix< Type > > SuSp(const DimensionedField< scalar, volMesh > &, const GeometricField< Type, fvPatchField, volMesh > &)
Foam::RASModels::LaunderSharmaKE::f2
tmp< volScalarField > f2() const
Definition: LaunderSharmaKE.C:46
Foam::RASModels::LaunderSharmaKE::rhoField
BasicTurbulenceModel::rhoField rhoField
Definition: LaunderSharmaKE.H:128
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam
Namespace for OpenFOAM.
Definition: combustionModel.C:30
rho
rho
Definition: pEqn.H:3
Foam::RASModels::LaunderSharmaKE::kSource
virtual tmp< fvScalarMatrix > kSource() const
Definition: LaunderSharmaKE.C:65
Foam::RASModels::LaunderSharmaKE::LaunderSharmaKE
LaunderSharmaKE(const LaunderSharmaKE &)
Foam::solve
SolverPerformance< Type > solve(fvMatrix< Type > &, const dictionary &)
Solve returning the solution statistics given convergence tolerance.
Foam::fvm::ddt
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
divU
volScalarField divU(fvc::div(fvc::absolute(phi, U)))
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
Foam::fvc::magSqrGradGrad
tmp< volScalarField > magSqrGradGrad(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcMagSqrGradGrad.C:44
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::eddyViscosity
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Foam::RASModels::LaunderSharmaKE::correctNut
virtual void correctNut()
Definition: LaunderSharmaKE.C:55
Foam::type
fileName::Type type(const fileName &)
Return the file type: DIRECTORY or FILE.
Definition: POSIX.C:588
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::fvm::div
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:46
Foam::tmp::clear
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:172
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::twoSymm
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:93
Foam::magSqr
dimensioned< scalar > magSqr(const dimensioned< Type > &)
Foam::RASModels::LaunderSharmaKE::transportModel
BasicTurbulenceModel::transportModel transportModel
Definition: LaunderSharmaKE.H:129
Foam::fvc::absolute
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:187
Foam::dev
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
Definition: dimensionedSymmTensor.C:104