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-2021 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 "fvModels.H"
28 #include "fvConstraints.H"
29 #include "bound.H"
30 
31 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace RASModels
36 {
37 
38 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
39 
40 template<class BasicMomentumTransportModel>
42 {
43  return exp(-3.4/sqr(scalar(1) + sqr(k_)/(this->nu()*epsilon_)/50.0));
44 }
45 
46 
47 template<class BasicMomentumTransportModel>
49 {
50  return
51  scalar(1)
52  - 0.3*exp(-min(sqr(sqr(k_)/(this->nu()*epsilon_)), scalar(50.0)));
53 }
54 
55 
56 template<class BasicMomentumTransportModel>
58 {
59  this->nut_ = Cmu_*fMu()*sqr(k_)/epsilon_;
60  this->nut_.correctBoundaryConditions();
61  fvConstraints::New(this->mesh_).constrain(this->nut_);
62 }
63 
64 
65 template<class BasicMomentumTransportModel>
68 {
69  return tmp<fvScalarMatrix>
70  (
71  new fvScalarMatrix
72  (
73  k_,
74  dimVolume*this->rho_.dimensions()*k_.dimensions()
75  /dimTime
76  )
77  );
78 }
79 
80 
81 template<class BasicMomentumTransportModel>
84 {
85  return tmp<fvScalarMatrix>
86  (
87  new fvScalarMatrix
88  (
89  epsilon_,
90  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
91  /dimTime
92  )
93  );
94 }
95 
96 
97 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
98 
99 template<class BasicMomentumTransportModel>
101 (
102  const alphaField& alpha,
103  const rhoField& rho,
104  const volVectorField& U,
105  const surfaceScalarField& alphaRhoPhi,
106  const surfaceScalarField& phi,
107  const viscosity& viscosity,
108  const word& type
109 )
110 :
112  (
113  type,
114  alpha,
115  rho,
116  U,
117  alphaRhoPhi,
118  phi,
119  viscosity
120  ),
121 
122  Cmu_
123  (
125  (
126  "Cmu",
127  this->coeffDict_,
128  0.09
129  )
130  ),
131  C1_
132  (
134  (
135  "C1",
136  this->coeffDict_,
137  1.44
138  )
139  ),
140  C2_
141  (
143  (
144  "C2",
145  this->coeffDict_,
146  1.92
147  )
148  ),
149  C3_
150  (
152  (
153  "C3",
154  this->coeffDict_,
155  0
156  )
157  ),
158  sigmak_
159  (
161  (
162  "sigmak",
163  this->coeffDict_,
164  1.0
165  )
166  ),
167  sigmaEps_
168  (
170  (
171  "sigmaEps",
172  this->coeffDict_,
173  1.3
174  )
175  ),
176 
177  k_
178  (
179  IOobject
180  (
181  IOobject::groupName("k", alphaRhoPhi.group()),
182  this->runTime_.timeName(),
183  this->mesh_,
186  ),
187  this->mesh_
188  ),
189 
190  epsilon_
191  (
192  IOobject
193  (
194  IOobject::groupName("epsilon", alphaRhoPhi.group()),
195  this->runTime_.timeName(),
196  this->mesh_,
199  ),
200  this->mesh_
201  )
202 {
203  bound(k_, this->kMin_);
204  bound(epsilon_, this->epsilonMin_);
205 
206  if (type == typeName)
207  {
208  this->printCoeffs(type);
209  }
210 }
211 
212 
213 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
214 
215 template<class BasicMomentumTransportModel>
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 BasicMomentumTransportModel>
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  const Foam::fvModels& fvModels(Foam::fvModels::New(this->mesh_));
252  (
253  Foam::fvConstraints::New(this->mesh_)
254  );
255 
257 
259 
260  // Calculate parameters and coefficients for Launder-Sharma low-Reynolds
261  // number model
262 
263  volScalarField E(2.0*this->nu()*nut*fvc::magSqrGradGrad(U));
264  volScalarField D(2.0*this->nu()*magSqr(fvc::grad(sqrt(k_))));
265 
266  tmp<volTensorField> tgradU = fvc::grad(U);
267  volScalarField G(this->GName(), nut*(tgradU() && dev(twoSymm(tgradU()))));
268  tgradU.clear();
269 
270 
271  // Dissipation equation
272  tmp<fvScalarMatrix> epsEqn
273  (
274  fvm::ddt(alpha, rho, epsilon_)
275  + fvm::div(alphaRhoPhi, epsilon_)
276  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
277  ==
278  C1_*alpha*rho*G*epsilon_/k_
279  - fvm::SuSp(((2.0/3.0)*C1_ - C3_)*alpha*rho*divU, epsilon_)
280  - fvm::Sp(C2_*f2()*alpha*rho*epsilon_/k_, epsilon_)
281  + alpha*rho*E
282  + epsilonSource()
283  + fvModels.source(alpha, rho, epsilon_)
284  );
285 
286  epsEqn.ref().relax();
287  fvConstraints.constrain(epsEqn.ref());
288  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
289  solve(epsEqn);
290  fvConstraints.constrain(epsilon_);
291  bound(epsilon_, this->epsilonMin_);
292 
293 
294  // Turbulent kinetic energy equation
296  (
297  fvm::ddt(alpha, rho, k_)
298  + fvm::div(alphaRhoPhi, k_)
299  - fvm::laplacian(alpha*rho*DkEff(), k_)
300  ==
301  alpha*rho*G - fvm::SuSp(2.0/3.0*alpha*rho*divU, k_)
302  - fvm::Sp(alpha*rho*(epsilon_ + D)/k_, k_)
303  + kSource()
304  + fvModels.source(alpha, rho, k_)
305  );
306 
307  kEqn.ref().relax();
308  fvConstraints.constrain(kEqn.ref());
309  solve(kEqn);
310  fvConstraints.constrain(k_);
311  bound(k_, this->kMin_);
312 
313  correctNut();
314 }
315 
316 
317 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
318 
319 } // End namespace RASModels
320 } // End namespace Foam
321 
322 // ************************************************************************* //
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:237
U
Definition: pEqn.H:72
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/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
virtual bool read()
Re-read model coefficients if they have changed.
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:47
volScalarField alpha(IOobject("alpha", runTime.timeName(), mesh, IOobject::READ_IF_PRESENT, IOobject::AUTO_WRITE), lambda *max(Ua &U, zeroSensitivity))
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.
const dimensionedScalar G
Newtonian constant of gravitation.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const GeometricField< Type, fvPatchField, volMesh > &)
Eddy viscosity turbulence model base class.
Definition: eddyViscosity.H:49
virtual tmp< fvScalarMatrix > kSource() const
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:53
BasicMomentumTransportModel::rhoField rhoField
const dimensionSet dimTime
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)
Foam::fvConstraints & fvConstraints
A class for handling words, derived from string.
Definition: word.H:59
static word groupName(Name name, const word &group)
Finite volume models.
Definition: fvModels.H:60
Info<< "Reading field p\"<< endl;volScalarField p(IOobject("p", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);pressureReference pressureReference(p, simple.dict());mesh.schemes().setFluxRequired(p.name());Info<< "Reading field pa\"<< endl;volScalarField pa(IOobject("pa", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);Info<< "Reading field Ua\"<< endl;volVectorField Ua(IOobject("Ua", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);# 65 "/home/ubuntu/OpenFOAM-10/applications/solvers/incompressible/adjointShapeOptimisationFoam/createFields.H" 2label paRefCell=0;scalar paRefValue=0.0;setRefCell(pa, simple.dict(), paRefCell, paRefValue);mesh.schemes().setFluxRequired(pa.name());autoPtr< viscosityModel > viscosity(viscosityModel::New(mesh))
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:46
tmp< volScalarField > f2() const
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
phi
Definition: correctPhi.H:3
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.
static autoPtr< dictionary > New(Istream &)
Construct top-level dictionary on freestore from Istream.
Definition: dictionaryIO.C:96
tmp< volScalarField > divU
Definition: alphaSuSp.H:9
bool constrain(fvMatrix< Type > &eqn) const
Apply constraints to an equation.
Foam::fvModels & fvModels
tmp< surfaceScalarField > absolute(const tmp< surfaceScalarField > &tphi, const volVectorField &U)
Return the given relative flux in absolute form.
Definition: fvcMeshPhi.C:202
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.
volScalarField & bound(volScalarField &, const dimensionedScalar &lowerBound)
Bound the given scalar field if it has gone unbounded.
Definition: bound.C:33
Finite volume constraints.
Definition: fvConstraints.H:57
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 dimVolume
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
tmp< volScalarField > fMu() const
A class for managing temporary objects.
Definition: PtrList.H:53
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:98
scalar nut
Namespace for OpenFOAM.