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-2020 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 BasicMomentumTransportModel>
41 {
42  return exp(-3.4/sqr(scalar(1) + sqr(k_)/(this->nu()*epsilon_)/50.0));
43 }
44 
45 
46 template<class BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
57 {
58  this->nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
59  this->nut_.correctBoundaryConditions();
60  fv::options::New(this->mesh_).correct(this->nut_);
61 }
62 
63 
64 template<class BasicMomentumTransportModel>
67 {
68  return tmp<fvScalarMatrix>
69  (
70  new fvScalarMatrix
71  (
72  k_,
73  dimVolume*this->rho_.dimensions()*k_.dimensions()
74  /dimTime
75  )
76  );
77 }
78 
79 
80 template<class BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
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& type
108 )
109 :
111  (
112  type,
113  alpha,
114  rho,
115  U,
116  alphaRhoPhi,
117  phi,
118  transport
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
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 BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
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  fv::options& fvOptions(fv::options::New(this->mesh_));
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  + fvOptions(alpha, rho, epsilon_)
279  );
280 
281  epsEqn.ref().relax();
282  fvOptions.constrain(epsEqn.ref());
283  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
284  solve(epsEqn);
285  fvOptions.correct(epsilon_);
286  bound(epsilon_, this->epsilonMin_);
287 
288 
289  // Turbulent kinetic energy equation
291  (
292  fvm::ddt(alpha, rho, k_)
293  + fvm::div(alphaRhoPhi, k_)
294  - fvm::laplacian(alpha*rho*DkEff(), k_)
295  ==
296  alpha*rho*G - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
297  - fvm::Sp(alpha*rho*(epsilon_ + D)/k_, k_)
298  + kSource()
299  + fvOptions(alpha, rho, k_)
300  );
301 
302  kEqn.ref().relax();
303  fvOptions.constrain(kEqn.ref());
304  solve(kEqn);
305  fvOptions.correct(k_);
306  bound(k_, this->kMin_);
307 
308  correctNut();
309 }
310 
311 
312 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
313 
314 } // End namespace RASModels
315 } // End namespace Foam
316 
317 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
BasicMomentumTransportModel::alphaField alphaField
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
virtual bool read()
Re-read model coefficients if they have changed.
void correct(GeometricField< Type, fvPatchField, volMesh > &field)
Apply correction to field.
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)
phi
Definition: pEqn.H:104
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
BasicMomentumTransportModel::rhoField rhoField
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
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 > &)
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(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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:68
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
LaunderSharmaKE(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const transportModel &transport, const word &type=typeName)
Construct from components.
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
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
BasicMomentumTransportModel::transportModel transportModel
const dimensionedScalar & G
Newtonian constant of gravitation.
scalar nut
Namespace for OpenFOAM.