kineticTheoryModel.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019-2023 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "kineticTheoryModel.H"
30 #include "mathematicalConstants.H"
31 #include "twoPhaseSystem.H"
32 #include "fvOptions.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 Foam::RASModels::kineticTheoryModel::kineticTheoryModel
37 (
38  const volScalarField& alpha,
39  const volScalarField& rho,
40  const volVectorField& U,
41  const surfaceScalarField& alphaRhoPhi,
42  const surfaceScalarField& phi,
43  const transportModel& phase,
44  const word& propertiesName,
45  const word& type
46 )
47 :
48  eddyViscosity
49  <
51  >
52  (
53  type,
54  alpha,
55  rho,
56  U,
57  alphaRhoPhi,
58  phi,
59  phase,
60  propertiesName
61  ),
62 
63  phase_(phase),
64 
65  viscosityModel_
66  (
67  kineticTheoryModels::viscosityModel::New
68  (
69  coeffDict_
70  )
71  ),
72  conductivityModel_
73  (
74  kineticTheoryModels::conductivityModel::New
75  (
76  coeffDict_
77  )
78  ),
79  radialModel_
80  (
81  kineticTheoryModels::radialModel::New
82  (
83  coeffDict_
84  )
85  ),
86  granularPressureModel_
87  (
88  kineticTheoryModels::granularPressureModel::New
89  (
90  coeffDict_
91  )
92  ),
93  frictionalStressModel_
94  (
95  kineticTheoryModels::frictionalStressModel::New
96  (
97  coeffDict_
98  )
99  ),
100 
101  equilibrium_(coeffDict_.get<bool>("equilibrium")),
102  e_("e", dimless, coeffDict_),
103  alphaMax_("alphaMax", dimless, coeffDict_),
104  alphaMinFriction_("alphaMinFriction", dimless, coeffDict_),
105  residualAlpha_("residualAlpha", dimless, coeffDict_),
106  maxNut_("maxNut", dimViscosity, 1000, coeffDict_),
107 
108  Theta_
109  (
110  IOobject
111  (
112  IOobject::groupName("Theta", phase.name()),
113  U.time().timeName(),
114  U.mesh(),
115  IOobject::MUST_READ,
116  IOobject::AUTO_WRITE
117  ),
118  U.mesh()
119  ),
120 
121  lambda_
122  (
123  IOobject
124  (
125  IOobject::groupName("lambda", phase.name()),
126  U.time().timeName(),
127  U.mesh(),
128  IOobject::NO_READ,
129  IOobject::NO_WRITE
130  ),
131  U.mesh(),
133  ),
134 
135  gs0_
136  (
137  IOobject
138  (
139  IOobject::groupName("gs0", phase.name()),
140  U.time().timeName(),
141  U.mesh(),
142  IOobject::NO_READ,
143  IOobject::NO_WRITE
144  ),
145  U.mesh(),
147  ),
148 
149  kappa_
150  (
151  IOobject
152  (
153  IOobject::groupName("kappa", phase.name()),
154  U.time().timeName(),
155  U.mesh(),
156  IOobject::NO_READ,
157  IOobject::NO_WRITE
158  ),
159  U.mesh(),
160  dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
161  ),
162 
163  nuFric_
164  (
165  IOobject
166  (
167  IOobject::groupName("nuFric", phase.name()),
168  U.time().timeName(),
169  U.mesh(),
170  IOobject::NO_READ,
171  IOobject::AUTO_WRITE
172  ),
173  U.mesh(),
175  )
176 {
177  if (type == typeName)
178  {
179  printCoeffs(type);
180  }
181 }
182 
183 
184 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
185 
187 {}
188 
189 
190 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
191 
193 {
194  if
195  (
196  eddyViscosity
197  <
198  RASModel<EddyDiffusivity<phaseCompressibleTurbulenceModel>>
199  >::read()
200  )
201  {
202  coeffDict().readEntry("equilibrium", equilibrium_);
203  e_.readIfPresent(coeffDict());
204  alphaMax_.readIfPresent(coeffDict());
205  alphaMinFriction_.readIfPresent(coeffDict());
206 
207  viscosityModel_->read();
208  conductivityModel_->read();
209  radialModel_->read();
210  granularPressureModel_->read();
211  frictionalStressModel_->read();
212 
213  return true;
214  }
215 
216  return false;
217 }
218 
219 
222 {
224  return nut_;
225 }
226 
227 
230 {
232  return nut_;
233 }
234 
235 
238 {
240  return nullptr;
241 }
242 
243 
246 {
247  return tmp<volSymmTensorField>
248  (
250  (
251  IOobject
252  (
253  IOobject::groupName("R", U_.group()),
254  runTime_.timeName(),
255  mesh_,
258  ),
259  - (nut_)*devTwoSymm(fvc::grad(U_))
260  - (lambda_*fvc::div(phi_))*symmTensor::I
261  )
262  );
263 }
264 
265 
268 {
269  const volScalarField& rho = phase_.rho();
270 
271  tmp<volScalarField> tpPrime
272  (
273  Theta_
274  *granularPressureModel_->granularPressureCoeffPrime
275  (
276  alpha_,
277  radialModel_->g0(alpha_, alphaMinFriction_, alphaMax_),
278  radialModel_->g0prime(alpha_, alphaMinFriction_, alphaMax_),
279  rho,
280  e_
281  )
282  + frictionalStressModel_->frictionalPressurePrime
283  (
284  phase_,
285  alphaMinFriction_,
286  alphaMax_
287  )
288  );
289 
290  volScalarField::Boundary& bpPrime =
291  tpPrime.ref().boundaryFieldRef();
292 
293  forAll(bpPrime, patchi)
294  {
295  if (!bpPrime[patchi].coupled())
296  {
297  bpPrime[patchi] == 0;
298  }
299  }
300 
301  return tpPrime;
302 }
303 
304 
307 {
308  return fvc::interpolate(pPrime());
309 }
310 
311 
314 {
315  return devRhoReff(U_);
316 }
317 
318 
321 (
322  const volVectorField& U
323 ) const
324 {
325  return tmp<volSymmTensorField>
326  (
328  (
329  IOobject
330  (
331  IOobject::groupName("devRhoReff", U.group()),
332  runTime_.timeName(),
333  mesh_,
336  ),
337  - (rho_*nut_)
339  - ((rho_*lambda_)*fvc::div(phi_))*symmTensor::I
340  )
341  );
342 }
343 
344 
347 (
349 ) const
350 {
351  return
352  (
353  - fvm::laplacian(rho_*nut_, U)
354  - fvc::div
355  (
356  (rho_*nut_)*dev2(T(fvc::grad(U)))
357  + ((rho_*lambda_)*fvc::div(phi_))
358  *dimensioned<symmTensor>("I", dimless, symmTensor::I)
359  )
360  );
361 }
362 
363 
365 {
366  // Local references
367  const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(phase_.fluid());
368  volScalarField alpha(max(alpha_, scalar(0)));
369  const volScalarField& rho = phase_.rho();
370  const surfaceScalarField& alphaRhoPhi = alphaRhoPhi_;
371  const volVectorField& U = U_;
372  const volVectorField& Uc_ = fluid.otherPhase(phase_).U();
373 
374  const scalar sqrtPi = sqrt(constant::mathematical::pi);
375  dimensionedScalar ThetaSmall("ThetaSmall", Theta_.dimensions(), 1.0e-6);
376  dimensionedScalar ThetaSmallSqrt(sqrt(ThetaSmall));
377 
378  tmp<volScalarField> tda(phase_.d());
379  const volScalarField& da = tda();
380 
381  tmp<volTensorField> tgradU(fvc::grad(U_));
382  const volTensorField& gradU(tgradU());
383  volSymmTensorField D(symm(gradU));
384 
385  // Calculating the radial distribution function
386  gs0_ = radialModel_->g0(alpha, alphaMinFriction_, alphaMax_);
387 
388  if (!equilibrium_)
389  {
390  // Particle viscosity (Table 3.2, p.47)
391  nut_ = viscosityModel_->nu(alpha, Theta_, gs0_, rho, da, e_);
392 
393  volScalarField ThetaSqrt("sqrtTheta", sqrt(Theta_));
394 
395  // Bulk viscosity p. 45 (Lun et al. 1984).
396  lambda_ = (4.0/3.0)*sqr(alpha)*da*gs0_*(1.0 + e_)*ThetaSqrt/sqrtPi;
397 
398  // Stress tensor, Definitions, Table 3.1, p. 43
400  (
401  rho*(2.0*nut_*D + (lambda_ - (2.0/3.0)*nut_)*tr(D)*I)
402  );
403 
404  // Dissipation (Eq. 3.24, p.50)
405  volScalarField gammaCoeff
406  (
407  "gammaCoeff",
408  12.0*(1.0 - sqr(e_))
409  *max(sqr(alpha), residualAlpha_)
410  *rho*gs0_*(1.0/da)*ThetaSqrt/sqrtPi
411  );
412 
413  // Drag
415  (
416  fluid.lookupSubModel<dragModel>
417  (
418  phase_,
419  fluid.otherPhase(phase_)
420  ).K()
421  );
422 
423  // Eq. 3.25, p. 50 Js = J1 - J2
424  volScalarField J1("J1", 3.0*beta);
425  volScalarField J2
426  (
427  "J2",
428  0.25*sqr(beta)*da*magSqr(U - Uc_)
429  /(
430  max(alpha, residualAlpha_)*rho
431  *sqrtPi*(ThetaSqrt + ThetaSmallSqrt)
432  )
433  );
434 
435  // particle pressure - coefficient in front of Theta (Eq. 3.22, p. 45)
436  volScalarField PsCoeff
437  (
438  granularPressureModel_->granularPressureCoeff
439  (
440  alpha,
441  gs0_,
442  rho,
443  e_
444  )
445  );
446 
447  // 'thermal' conductivity (Table 3.3, p. 49)
448  kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rho, da, e_);
449 
450  fv::options& fvOptions(fv::options::New(mesh_));
451 
452  // Construct the granular temperature equation (Eq. 3.20, p. 44)
453  // NB. note that there are two typos in Eq. 3.20:
454  // Ps should be without grad
455  // the laplacian has the wrong sign
456  fvScalarMatrix ThetaEqn
457  (
458  1.5*
459  (
460  fvm::ddt(alpha, rho, Theta_)
461  + fvm::div(alphaRhoPhi, Theta_)
462  - fvc::Sp(fvc::ddt(alpha, rho) + fvc::div(alphaRhoPhi), Theta_)
463  )
464  - fvm::laplacian(kappa_, Theta_, "laplacian(kappa,Theta)")
465  ==
466  - fvm::SuSp((PsCoeff*I) && gradU, Theta_)
467  + (tau && gradU)
468  + fvm::Sp(-gammaCoeff, Theta_)
469  + fvm::Sp(-J1, Theta_)
470  + fvm::Sp(J2/(Theta_ + ThetaSmall), Theta_)
471  + fvOptions(alpha, rho, Theta_)
472  );
473 
474  ThetaEqn.relax();
475  fvOptions.constrain(ThetaEqn);
476  ThetaEqn.solve();
477  fvOptions.correct(Theta_);
478  }
479  else
480  {
481  // Equilibrium => dissipation == production
482  // Eq. 4.14, p.82
483  volScalarField K1("K1", 2.0*(1.0 + e_)*rho*gs0_);
485  (
486  "K3",
487  0.5*da*rho*
488  (
489  (sqrtPi/(3.0*(3.0 - e_)))
490  *(1.0 + 0.4*(1.0 + e_)*(3.0*e_ - 1.0)*alpha*gs0_)
491  +1.6*alpha*gs0_*(1.0 + e_)/sqrtPi
492  )
493  );
494 
496  (
497  "K2",
498  4.0*da*rho*(1.0 + e_)*alpha*gs0_/(3.0*sqrtPi) - 2.0*K3/3.0
499  );
500 
501  volScalarField K4("K4", 12.0*(1.0 - sqr(e_))*rho*gs0_/(da*sqrtPi));
502 
503  volScalarField trD
504  (
505  "trD",
506  alpha/(alpha + residualAlpha_)
507  *fvc::div(phi_)
508  );
509  volScalarField tr2D("tr2D", sqr(trD));
510  volScalarField trD2("trD2", tr(D & D));
511 
512  volScalarField t1("t1", K1*alpha + rho);
513  volScalarField l1("l1", -t1*trD);
514  volScalarField l2("l2", sqr(t1)*tr2D);
515  volScalarField l3
516  (
517  "l3",
518  4.0
519  *K4
520  *alpha
521  *(2.0*K3*trD2 + K2*tr2D)
522  );
523 
524  Theta_ = sqr
525  (
526  (l1 + sqrt(l2 + l3))
527  /(2.0*max(alpha, residualAlpha_)*K4)
528  );
529 
530  kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rho, da, e_);
531  }
532 
533  Theta_.clamp_range(0, 100);
534 
535  {
536  // particle viscosity (Table 3.2, p.47)
537  nut_ = viscosityModel_->nu(alpha, Theta_, gs0_, rho, da, e_);
538 
539  volScalarField ThetaSqrt("sqrtTheta", sqrt(Theta_));
540 
541  // Bulk viscosity p. 45 (Lun et al. 1984).
542  lambda_ = (4.0/3.0)*sqr(alpha)*da*gs0_*(1.0 + e_)*ThetaSqrt/sqrtPi;
543 
544  // Frictional pressure
545  volScalarField pf
546  (
547  frictionalStressModel_->frictionalPressure
548  (
549  phase_,
550  alphaMinFriction_,
551  alphaMax_
552  )
553  );
554 
555  nuFric_ = frictionalStressModel_->nu
556  (
557  phase_,
558  alphaMinFriction_,
559  alphaMax_,
560  pf/rho,
561  D
562  );
563 
564  // Limit viscosity and add frictional viscosity
565  nut_.clamp_max(maxNut_);
566 
567  nuFric_ = min(nuFric_, maxNut_ - nut_);
568  nut_ += nuFric_;
569  }
570 
571  if (debug)
572  {
573  Info<< typeName << ':' << nl
574  << " max(Theta) = " << max(Theta_).value() << nl
575  << " max(nut) = " << max(nut_).value() << endl;
576  }
577 }
578 
579 
580 // ************************************************************************* //
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
twoPhaseSystem & fluid
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
Definition: volFieldsFwd.H:85
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:37
type
Types of root.
Definition: Roots.H:52
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
#define K1
Definition: SHA1.C:143
#define K4
Definition: SHA1.C:146
#define K2
Definition: SHA1.C:144
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
GeometricField< tensor, fvPatchField, volMesh > volTensorField
Definition: volFieldsFwd.H:86
dimensionedSymmTensor sqr(const dimensionedVector &dv)
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
constexpr char nl
The newline &#39;\n&#39; character (0x0a)
Definition: Ostream.H:49
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
#define K3
Definition: SHA1.C:145
bool coupled(solutionDict.getOrDefault("coupledEnergyField", false))
const dimensionSet dimViscosity
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:487
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tf1, const word &name, const dimensionSet &dimensions, const bool initCopy=false)
Global function forwards to reuseTmpDimensionedField::New.
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:62
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
Type of boundary fields.
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
GeometricField< vector, fvPatchField, volMesh > volVectorField
Definition: volFieldsFwd.H:82
CGAL::Exact_predicates_exact_constructions_kernel K
virtual bool read()
Re-read model coefficients if they have changed.
fv::options & fvOptions
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition: int32.H:127
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:414
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition: fvcDdt.C:40
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
word timeName
Definition: getTimeIndex.H:3
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition: POSIX.C:799
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
dynamicFvMesh & mesh
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:52
static const SymmTensor I
Definition: SymmTensor.H:74
static const Identity< scalar > I
Definition: Identity.H:100
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure&#39;.
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvmDdt.C:41
constexpr scalar pi(M_PI)
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
virtual void printCoeffs(const word &type)
Print model coefficients.
Definition: RASModel.C:27
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition: fvmDiv.C:41
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
int debug
Static debugging option.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
static tmp< GeometricField< Type, fvsPatchField, surfaceMesh > > interpolate(const GeometricField< Type, fvPatchField, volMesh > &tvf, const surfaceScalarField &faceFlux, Istream &schemeData)
Interpolate field onto faces using scheme given by Istream.
void clamp_max(const Type &upper)
Impose upper (ceiling) clamp on the field values (in-place)
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure&#39;.
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
U
Definition: pEqn.H:72
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:41
Nothing to be read.
messageStream Info
Information stream (stdout output on master, null elsewhere)
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
const dimensionedScalar & D
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present.
Definition: fvOptions.C:96
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition: error.H:666
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
Definition: SymmTensorI.H:491
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:133