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-2023 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 "fvcMagSqrGradGrad.H"
28 #include "fvModels.H"
29 #include "fvConstraints.H"
30 #include "bound.H"
31 
32 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace RASModels
37 {
38 
39 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
40 
41 template<class BasicMomentumTransportModel>
43 {
44  return exp(-3.4/sqr(scalar(1) + sqr(k_)/(this->nu()*epsilon_)/50.0));
45 }
46 
47 
48 template<class BasicMomentumTransportModel>
50 {
51  return
52  scalar(1)
53  - 0.3*exp(-min(sqr(sqr(k_)/(this->nu()*epsilon_)), scalar(50.0)));
54 }
55 
56 
57 template<class BasicMomentumTransportModel>
59 {
60  this->nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
61  this->nut_.correctBoundaryConditions();
62  fvConstraints::New(this->mesh_).constrain(this->nut_);
63 }
64 
65 
66 template<class BasicMomentumTransportModel>
69 {
70  return tmp<fvScalarMatrix>
71  (
72  new fvScalarMatrix
73  (
74  k_,
75  dimVolume*this->rho_.dimensions()*k_.dimensions()
76  /dimTime
77  )
78  );
79 }
80 
81 
82 template<class BasicMomentumTransportModel>
85 {
86  return tmp<fvScalarMatrix>
87  (
88  new fvScalarMatrix
89  (
90  epsilon_,
91  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
92  /dimTime
93  )
94  );
95 }
96 
97 
98 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
99 
100 template<class BasicMomentumTransportModel>
102 (
103  const alphaField& alpha,
104  const rhoField& rho,
105  const volVectorField& U,
106  const surfaceScalarField& alphaRhoPhi,
107  const surfaceScalarField& phi,
108  const viscosity& viscosity,
109  const word& type
110 )
111 :
112  eddyViscosity<RASModel<BasicMomentumTransportModel>>
113  (
114  type,
115  alpha,
116  rho,
117  U,
118  alphaRhoPhi,
119  phi,
120  viscosity
121  ),
122 
123  Cmu_
124  (
125  dimensioned<scalar>::lookupOrAddToDict
126  (
127  "Cmu",
128  this->coeffDict_,
129  0.09
130  )
131  ),
132  C1_
133  (
134  dimensioned<scalar>::lookupOrAddToDict
135  (
136  "C1",
137  this->coeffDict_,
138  1.44
139  )
140  ),
141  C2_
142  (
143  dimensioned<scalar>::lookupOrAddToDict
144  (
145  "C2",
146  this->coeffDict_,
147  1.92
148  )
149  ),
150  C3_
151  (
152  dimensioned<scalar>::lookupOrAddToDict
153  (
154  "C3",
155  this->coeffDict_,
156  0
157  )
158  ),
159  sigmak_
160  (
161  dimensioned<scalar>::lookupOrAddToDict
162  (
163  "sigmak",
164  this->coeffDict_,
165  1.0
166  )
167  ),
168  sigmaEps_
169  (
170  dimensioned<scalar>::lookupOrAddToDict
171  (
172  "sigmaEps",
173  this->coeffDict_,
174  1.3
175  )
176  ),
177 
178  k_
179  (
180  IOobject
181  (
182  this->groupName("k"),
183  this->runTime_.name(),
184  this->mesh_,
185  IOobject::MUST_READ,
186  IOobject::AUTO_WRITE
187  ),
188  this->mesh_
189  ),
190 
191  epsilon_
192  (
193  IOobject
194  (
195  this->groupName("epsilon"),
196  this->runTime_.name(),
197  this->mesh_,
198  IOobject::MUST_READ,
199  IOobject::AUTO_WRITE
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 BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
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  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
253  (
254  Foam::fvConstraints::New(this->mesh_)
255  );
256 
258 
259  volScalarField divU(fvc::div(fvc::absolute(this->phi(), U)));
260 
261  // Calculate parameters and coefficients for Launder-Sharma low-Reynolds
262  // number model
263 
264  volScalarField E(2.0*this->nu()*nut*fvc::magSqrGradGrad(U));
265  volScalarField D(2.0*this->nu()*magSqr(fvc::grad(sqrt(k_))));
266 
267  tmp<volTensorField> tgradU = fvc::grad(U);
268  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
269  tgradU.clear();
270 
271 
272  // Dissipation equation
273  tmp<fvScalarMatrix> epsEqn
274  (
275  fvm::ddt(alpha, rho, epsilon_)
276  + fvm::div(alphaRhoPhi, epsilon_)
277  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
278  ==
279  C1_*alpha*rho*G*epsilon_/k_
280  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha*rho*divU, epsilon_)
281  - fvm::Sp(C2_*f2()*alpha*rho*epsilon_/k_, epsilon_)
282  + alpha*rho*E
283  + epsilonSource()
284  + fvModels.source(alpha, rho, epsilon_)
285  );
286 
287  epsEqn.ref().relax();
288  fvConstraints.constrain(epsEqn.ref());
289  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
290  solve(epsEqn);
291  fvConstraints.constrain(epsilon_);
292  bound(epsilon_, this->epsilonMin_);
293 
294 
295  // Turbulent kinetic energy equation
297  (
298  fvm::ddt(alpha, rho, k_)
299  + fvm::div(alphaRhoPhi, k_)
300  - fvm::laplacian(alpha*rho*DkEff(), k_)
301  ==
302  alpha*rho*G - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
303  - fvm::Sp(alpha*rho*(epsilon_ + D)/k_, k_)
304  + kSource()
305  + fvModels.source(alpha, rho, k_)
306  );
307 
308  kEqn.ref().relax();
309  fvConstraints.constrain(kEqn.ref());
310  solve(kEqn);
312  bound(k_, this->kMin_);
313 
314  correctNut();
315 }
316 
317 
318 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
319 
320 } // End namespace RASModels
321 } // End namespace Foam
322 
323 // ************************************************************************* //
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
Bound the given scalar field where it is below the specified minimum.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:56
virtual tmp< fvScalarMatrix > epsilonSource() const
tmp< volScalarField > f2() const
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
LaunderSharmaKE(const alphaField &alpha, const rhoField &rho, const volVectorField &U, const surfaceScalarField &alphaRhoPhi, const surfaceScalarField &phi, const viscosity &viscosity, const word &type=typeName)
Construct from components.
virtual tmp< fvScalarMatrix > kSource() const
tmp< volScalarField > fMu() const
virtual bool read()
Re-read model coefficients if they have changed.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:100
Generic dimensioned Type class.
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite volume constraints.
Definition: fvConstraints.H:62
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
Finite volume models.
Definition: fvModels.H:65
BasicMomentumTransportModel::alphaField alphaField
BasicMomentumTransportModel::rhoField rhoField
A class for managing temporary objects.
Definition: tmp.H:55
void clear() const
If object pointer points to valid object:
Definition: tmpI.H:237
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Abstract base class for all fluid physical properties.
Definition: viscosity.H:50
A class for handling words, derived from string.
Definition: word.H:62
Foam::fvConstraints & fvConstraints(Foam::fvConstraints::New(mesh))
Foam::fvModels & fvModels(Foam::fvModels::New(mesh))
const scalar nut
Calculate the magnitude of the square of the gradient of the gradient of the given volField.
U
Definition: pEqn.H:72
volScalarField alpha(IOobject("alpha", runTime.name(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< volScalarField > magSqrGradGrad(const VolField< Type > &vf)
tmp< VolField< typename outerProduct< vector, Type >::type > > grad(const SurfaceField< Type > &ssf)
Definition: fvcGrad.C:46
tmp< VolField< Type > > div(const SurfaceField< Type > &ssf)
Definition: fvcDiv.C:47
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
tmp< fvMatrix< Type > > laplacian(const VolField< Type > &vf, const word &name)
Definition: fvmLaplacian.C:47
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
Namespace for OpenFOAM.
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimVolume
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensioned< scalar > magSqr(const dimensioned< Type > &)
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.