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-2024 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  epsilon_ = max(epsilon_, Cmu_*sqr(k_)/(this->nutMaxCoeff_*this->nu()));
61 }
62 
63 
64 template<class BasicMomentumTransportModel>
66 {
67  boundEpsilon();
68  this->nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
69  this->nut_.correctBoundaryConditions();
70  fvConstraints::New(this->mesh_).constrain(this->nut_);
71 }
72 
73 
74 template<class BasicMomentumTransportModel>
77 {
78  return tmp<fvScalarMatrix>
79  (
80  new fvScalarMatrix
81  (
82  k_,
83  dimVolume*this->rho_.dimensions()*k_.dimensions()
84  /dimTime
85  )
86  );
87 }
88 
89 
90 template<class BasicMomentumTransportModel>
93 {
94  return tmp<fvScalarMatrix>
95  (
96  new fvScalarMatrix
97  (
98  epsilon_,
99  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
100  /dimTime
101  )
102  );
103 }
104 
105 
106 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
107 
108 template<class BasicMomentumTransportModel>
110 (
111  const alphaField& alpha,
112  const rhoField& rho,
113  const volVectorField& U,
114  const surfaceScalarField& alphaRhoPhi,
115  const surfaceScalarField& phi,
116  const viscosity& viscosity,
117  const word& type
118 )
119 :
120  eddyViscosity<RASModel<BasicMomentumTransportModel>>
121  (
122  type,
123  alpha,
124  rho,
125  U,
126  alphaRhoPhi,
127  phi,
128  viscosity
129  ),
130 
131  Cmu_("Cmu", this->coeffDict(), 0.09),
132  C1_("C1", this->coeffDict(), 1.44),
133  C2_("C2", this->coeffDict(), 1.92),
134  C3_("C3", this->coeffDict(), 0),
135  sigmak_("sigmak", this->coeffDict(), 1.0),
136  sigmaEps_("sigmaEps", this->coeffDict(), 1.3),
137 
138  k_
139  (
140  IOobject
141  (
142  this->groupName("k"),
143  this->runTime_.name(),
144  this->mesh_,
145  IOobject::MUST_READ,
146  IOobject::AUTO_WRITE
147  ),
148  this->mesh_
149  ),
150 
151  epsilon_
152  (
153  IOobject
154  (
155  this->groupName("epsilon"),
156  this->runTime_.name(),
157  this->mesh_,
158  IOobject::MUST_READ,
159  IOobject::AUTO_WRITE
160  ),
161  this->mesh_
162  )
163 {
164  bound(k_, this->kMin_);
165  boundEpsilon();
166 }
167 
168 
169 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
170 
171 template<class BasicMomentumTransportModel>
173 {
175  {
176  Cmu_.readIfPresent(this->coeffDict());
177  C1_.readIfPresent(this->coeffDict());
178  C2_.readIfPresent(this->coeffDict());
179  C3_.readIfPresent(this->coeffDict());
180  sigmak_.readIfPresent(this->coeffDict());
181  sigmaEps_.readIfPresent(this->coeffDict());
182 
183  return true;
184  }
185  else
186  {
187  return false;
188  }
189 }
190 
191 
192 template<class BasicMomentumTransportModel>
194 {
195  if (!this->turbulence_)
196  {
197  return;
198  }
199 
200  // Local references
201  const alphaField& alpha = this->alpha_;
202  const rhoField& rho = this->rho_;
203  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
204  const volVectorField& U = this->U_;
205  volScalarField& nut = this->nut_;
206  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
208  (
209  Foam::fvConstraints::New(this->mesh_)
210  );
211 
213 
214  volScalarField divU(fvc::div(fvc::absolute(this->phi(), U)));
215 
216  // Calculate parameters and coefficients for Launder-Sharma low-Reynolds
217  // number model
218 
219  volScalarField E(2.0*this->nu()*nut*fvc::magSqrGradGrad(U));
220  volScalarField D(2.0*this->nu()*magSqr(fvc::grad(sqrt(k_))));
221 
222  tmp<volTensorField> tgradU = fvc::grad(U);
223  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
224  tgradU.clear();
225 
226 
227  // Dissipation equation
228  tmp<fvScalarMatrix> epsEqn
229  (
230  fvm::ddt(alpha, rho, epsilon_)
231  + fvm::div(alphaRhoPhi, epsilon_)
232  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
233  ==
234  C1_*alpha*rho*G*epsilon_/k_
235  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha*rho*divU, epsilon_)
236  - fvm::Sp(C2_*f2()*alpha*rho*epsilon_/k_, epsilon_)
237  + alpha*rho*E
238  + epsilonSource()
239  + fvModels.source(alpha, rho, epsilon_)
240  );
241 
242  epsEqn.ref().relax();
243  fvConstraints.constrain(epsEqn.ref());
244  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
245  solve(epsEqn);
246  fvConstraints.constrain(epsilon_);
247  boundEpsilon();
248 
249 
250  // Turbulent kinetic energy equation
252  (
253  fvm::ddt(alpha, rho, k_)
254  + fvm::div(alphaRhoPhi, k_)
255  - fvm::laplacian(alpha*rho*DkEff(), k_)
256  ==
257  alpha*rho*G - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
258  - fvm::Sp(alpha*rho*(epsilon_ + D)/k_, k_)
259  + kSource()
260  + fvModels.source(alpha, rho, k_)
261  );
262 
263  kEqn.ref().relax();
264  fvConstraints.constrain(kEqn.ref());
265  solve(kEqn);
267  bound(k_, this->kMin_);
268 
269  correctNut();
270 }
271 
272 
273 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
274 
275 } // End namespace RASModels
276 } // End namespace Foam
277 
278 // ************************************************************************* //
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
Source term for the epsilon equation.
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.
void boundEpsilon()
Bound epsilon.
virtual void correctNut()
Correct the eddy-viscosity nut.
virtual tmp< fvScalarMatrix > kSource() const
Source term for the k equation.
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:103
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:52
Finite volume constraints.
Definition: fvConstraints.H:67
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
tmp< fvMatrix< Type > > source(const VolField< Type > &field) const
Return source for an equation.
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:253
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:197
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 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 > > div(const surfaceScalarField &flux, const VolField< Type > &vf, const word &name)
Definition: fvmDiv.C:48
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > SuSp(const volScalarField::Internal &, const VolField< Type > &)
tmp< fvMatrix< Type > > ddt(const VolField< Type > &vf)
Definition: fvmDdt.C:46
static const coefficient D("D", dimTemperature, 257.14)
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
bool read(const char *, int32_t &)
Definition: int32IO.C:85
void dev(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
void twoSymm(LagrangianPatchField< tensor > &f, const LagrangianPatchField< tensor > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimTime
const dimensionSet dimVolume
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
bool bound(volScalarField &, const dimensionedScalar &min)
Bound the given scalar field where it is below the specified min value.
Definition: bound.C:31
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void magSqr(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
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.