realizableKE.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 "realizableKE.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  const volTensorField& gradU,
43  const volScalarField& S2,
44  const volScalarField& magS
45 )
46 {
47  tmp<volSymmTensorField> tS = dev(symm(gradU));
48  const volSymmTensorField& S = tS();
49 
51  (
52  (2*sqrt(2.0))*((S&S)&&S)
53  /(
54  magS*S2
55  + dimensionedScalar(dimensionSet(0, 0, -3, 0, 0), small)
56  )
57  );
58 
59  tS.clear();
60 
61  volScalarField phis
62  (
63  (1.0/3.0)*acos(min(max(sqrt(6.0)*W, -scalar(1)), scalar(1)))
64  );
65  volScalarField As(sqrt(6.0)*cos(phis));
66  volScalarField Us(sqrt(S2/2.0 + magSqr(skew(gradU))));
67 
68  return 1.0/(A0_ + As*Us*k_/epsilon_);
69 }
70 
71 
72 template<class BasicMomentumTransportModel>
74 (
75  const volTensorField& gradU,
76  const volScalarField& S2,
77  const volScalarField& magS
78 )
79 {
80  this->nut_ = rCmu(gradU, S2, magS)*sqr(k_)/epsilon_;
81  this->nut_.correctBoundaryConditions();
82  fv::options::New(this->mesh_).correct(this->nut_);
83 }
84 
85 
86 template<class BasicMomentumTransportModel>
88 {
89  tmp<volTensorField> tgradU = fvc::grad(this->U_);
90  volScalarField S2(modelName("S2"), 2*magSqr(dev(symm(tgradU()))));
91  volScalarField magS(modelName("magS"), sqrt(S2));
92 
93  correctNut(tgradU(), S2, magS);
94 }
95 
96 
97 template<class BasicMomentumTransportModel>
99 {
100  return tmp<fvScalarMatrix>
101  (
102  new fvScalarMatrix
103  (
104  k_,
105  dimVolume*this->rho_.dimensions()*k_.dimensions()
106  /dimTime
107  )
108  );
109 }
110 
111 
112 template<class BasicMomentumTransportModel>
115 {
116  return tmp<fvScalarMatrix>
117  (
118  new fvScalarMatrix
119  (
120  epsilon_,
121  dimVolume*this->rho_.dimensions()*epsilon_.dimensions()
122  /dimTime
123  )
124  );
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
129 
130 template<class BasicMomentumTransportModel>
132 (
133  const alphaField& alpha,
134  const rhoField& rho,
135  const volVectorField& U,
136  const surfaceScalarField& alphaRhoPhi,
137  const surfaceScalarField& phi,
138  const transportModel& transport,
139  const word& type
140 )
141 :
143  (
144  type,
145  alpha,
146  rho,
147  U,
148  alphaRhoPhi,
149  phi,
150  transport
151  ),
152  A0_
153  (
155  (
156  "A0",
157  this->coeffDict_,
158  4.0
159  )
160  ),
161  C2_
162  (
164  (
165  "C2",
166  this->coeffDict_,
167  1.9
168  )
169  ),
170  sigmak_
171  (
173  (
174  "sigmak",
175  this->coeffDict_,
176  1.0
177  )
178  ),
179  sigmaEps_
180  (
182  (
183  "sigmaEps",
184  this->coeffDict_,
185  1.2
186  )
187  ),
188 
189  k_
190  (
191  IOobject
192  (
193  "k",
194  this->runTime_.timeName(),
195  this->mesh_,
198  ),
199  this->mesh_
200  ),
201  epsilon_
202  (
203  IOobject
204  (
205  "epsilon",
206  this->runTime_.timeName(),
207  this->mesh_,
210  ),
211  this->mesh_
212  )
213 {
214  bound(k_, this->kMin_);
215  bound(epsilon_, this->epsilonMin_);
216 
217  if (type == typeName)
218  {
219  this->printCoeffs(type);
220  }
221 }
222 
223 
224 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
225 
226 template<class BasicMomentumTransportModel>
228 {
230  {
231  A0_.readIfPresent(this->coeffDict());
232  C2_.readIfPresent(this->coeffDict());
233  sigmak_.readIfPresent(this->coeffDict());
234  sigmaEps_.readIfPresent(this->coeffDict());
235 
236  return true;
237  }
238  else
239  {
240  return false;
241  }
242 }
243 
244 
245 template<class BasicMomentumTransportModel>
247 {
248  if (!this->turbulence_)
249  {
250  return;
251  }
252 
253  // Local references
254  const alphaField& alpha = this->alpha_;
255  const rhoField& rho = this->rho_;
256  const surfaceScalarField& alphaRhoPhi = this->alphaRhoPhi_;
257  const volVectorField& U = this->U_;
258  volScalarField& nut = this->nut_;
259  fv::options& fvOptions(fv::options::New(this->mesh_));
260 
262 
264  (
265  modelName("divU"),
266  fvc::div(fvc::absolute(this->phi(), U))()
267  );
268 
269  tmp<volTensorField> tgradU = fvc::grad(U);
270  volScalarField S2(modelName("S2"), 2*magSqr(dev(symm(tgradU()))));
271  volScalarField magS(modelName("magS"), sqrt(S2));
272 
273  volScalarField::Internal eta(modelName("eta"), magS()*k_()/epsilon_());
275  (
276  modelName("C1"),
277  max(eta/(scalar(5) + eta), scalar(0.43))
278  );
279 
281  (
282  this->GName(),
283  nut*(tgradU().v() && dev(twoSymm(tgradU().v())))
284  );
285 
286  // Update epsilon and G at the wall
287  epsilon_.boundaryFieldRef().updateCoeffs();
288 
289  // Dissipation equation
290  tmp<fvScalarMatrix> epsEqn
291  (
292  fvm::ddt(alpha, rho, epsilon_)
293  + fvm::div(alphaRhoPhi, epsilon_)
294  - fvm::laplacian(alpha*rho*DepsilonEff(), epsilon_)
295  ==
296  C1*alpha()*rho()*magS()*epsilon_()
297  - fvm::Sp
298  (
299  C2_*alpha()*rho()*epsilon_()/(k_() + sqrt(this->nu()()*epsilon_())),
300  epsilon_
301  )
302  + epsilonSource()
303  + fvOptions(alpha, rho, epsilon_)
304  );
305 
306  epsEqn.ref().relax();
307  fvOptions.constrain(epsEqn.ref());
308  epsEqn.ref().boundaryManipulate(epsilon_.boundaryFieldRef());
309  solve(epsEqn);
310  fvOptions.correct(epsilon_);
311  bound(epsilon_, this->epsilonMin_);
312 
313 
314  // Turbulent kinetic energy equation
315 
317  (
318  fvm::ddt(alpha, rho, k_)
319  + fvm::div(alphaRhoPhi, k_)
320  - fvm::laplacian(alpha*rho*DkEff(), k_)
321  ==
322  alpha()*rho()*G
323  - fvm::SuSp(2.0/3.0*alpha()*rho()*divU, k_)
324  - fvm::Sp(alpha()*rho()*epsilon_()/k_(), k_)
325  + kSource()
326  + fvOptions(alpha, rho, k_)
327  );
328 
329  kEqn.ref().relax();
330  fvOptions.constrain(kEqn.ref());
331  solve(kEqn);
332  fvOptions.correct(k_);
333  bound(k_, this->kMin_);
334 
335  correctNut(tgradU(), S2, magS);
336 }
337 
338 
339 // * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * * //
340 
341 } // End namespace RASModels
342 } // End namespace Foam
343 
344 // ************************************************************************* //
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
dimensionedScalar acos(const dimensionedScalar &ds)
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedTensor skew(const dimensionedTensor &dt)
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)
realizableKE(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.
Definition: realizableKE.C:132
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
virtual tmp< fvScalarMatrix > kSource() const
Definition: realizableKE.C:98
Finite-volume options.
Definition: fvOptions.H:52
BasicMomentumTransportModel::alphaField alphaField
Definition: RASModel.H:88
virtual bool read()
Re-read model coefficients if they have changed.
Definition: realizableKE.C:227
dimensionedSymmTensor twoSymm(const dimensionedSymmTensor &dt)
phi
Definition: pEqn.H:104
Templated abstract base class for RAS turbulence models.
Definition: RASModel.H:49
tmp< volScalarField > rCmu(const volTensorField &gradU, const volScalarField &S2, const volScalarField &magS)
Definition: realizableKE.C:41
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Dimension set for the base types.
Definition: dimensionSet.H:120
bool read(const char *, int32_t &)
Definition: int32IO.C:85
dimensionedSymmTensor dev(const dimensionedSymmTensor &dt)
dimensionedScalar cos(const dimensionedScalar &ds)
BasicMomentumTransportModel::transportModel transportModel
Definition: RASModel.H:90
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
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 > &)
virtual tmp< fvScalarMatrix > epsilonSource() const
Definition: realizableKE.C:114
BasicMomentumTransportModel::rhoField rhoField
Definition: RASModel.H:89
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
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmLaplacian.C:46
virtual void correct()
Solve the turbulence equations and correct the turbulence viscosity.
Definition: realizableKE.C:246
const scalarList W(::W(thermo))
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
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
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
const dimensionedScalar & G
Newtonian constant of gravitation.
scalar nut
Namespace for OpenFOAM.