Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
11License
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"
31#include "twoPhaseSystem.H"
32#include "fvOptions.H"
33
34// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35
36Foam::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 IOobject::REGISTER
118 ),
119 U.mesh()
120 ),
121
122 lambda_
123 (
124 IOobject
125 (
126 IOobject::groupName("lambda", phase.name()),
127 U.time().timeName(),
128 U.mesh(),
129 IOobject::NO_READ,
130 IOobject::NO_WRITE
131 ),
132 U.mesh(),
134 ),
135
136 gs0_
137 (
138 IOobject
139 (
140 IOobject::groupName("gs0", phase.name()),
141 U.time().timeName(),
142 U.mesh(),
143 IOobject::NO_READ,
144 IOobject::NO_WRITE
145 ),
146 U.mesh(),
148 ),
149
150 kappa_
151 (
152 IOobject
153 (
154 IOobject::groupName("kappa", phase.name()),
155 U.time().timeName(),
156 U.mesh(),
157 IOobject::NO_READ,
158 IOobject::NO_WRITE
159 ),
160 U.mesh(),
161 dimensionedScalar(dimensionSet(1, -1, -1, 0, 0), Zero)
162 ),
163
164 nuFric_
165 (
166 IOobject
167 (
168 IOobject::groupName("nuFric", phase.name()),
169 U.time().timeName(),
170 U.mesh(),
171 IOobject::NO_READ,
172 IOobject::AUTO_WRITE,
173 IOobject::REGISTER
174 ),
175 U.mesh(),
177 )
178{
179 if (type == typeName)
180 {
181 printCoeffs(type);
182 }
183}
184
185
186// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
187
189{}
190
191
192// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
193
195{
196 if
197 (
198 eddyViscosity
199 <
200 RASModel<EddyDiffusivity<phaseCompressibleTurbulenceModel>>
201 >::read()
202 )
203 {
204 coeffDict().readEntry("equilibrium", equilibrium_);
205 e_.readIfPresent(coeffDict());
206 alphaMax_.readIfPresent(coeffDict());
207 alphaMinFriction_.readIfPresent(coeffDict());
208
209 viscosityModel_->read();
210 conductivityModel_->read();
211 radialModel_->read();
212 granularPressureModel_->read();
213 frictionalStressModel_->read();
214
215 return true;
216 }
217
218 return false;
219}
220
221
222Foam::tmp<Foam::volScalarField>
224{
226 return nut_;
227}
228
229
230Foam::tmp<Foam::volScalarField>
232{
234 return nut_;
235}
236
237
238Foam::tmp<Foam::volScalarField>
240{
242 return nullptr;
243}
244
245
246Foam::tmp<Foam::volSymmTensorField>
248{
250 (
251 IOobject::groupName("R", U_.group()),
253 (
254 - (nut_)*devTwoSymm(fvc::grad(U_))
255 - (lambda_*fvc::div(phi_))*symmTensor::I
256 )
257 );
258}
259
260
261Foam::tmp<Foam::volScalarField>
263{
264 const volScalarField& rho = phase_.rho();
265
266 tmp<volScalarField> tpPrime
267 (
268 Theta_
269 *granularPressureModel_->granularPressureCoeffPrime
270 (
271 alpha_,
272 radialModel_->g0(alpha_, alphaMinFriction_, alphaMax_),
273 radialModel_->g0prime(alpha_, alphaMinFriction_, alphaMax_),
274 rho,
275 e_
276 )
277 + frictionalStressModel_->frictionalPressurePrime
278 (
279 phase_,
280 alphaMinFriction_,
281 alphaMax_
282 )
283 );
284
285 volScalarField::Boundary& bpPrime =
286 tpPrime.ref().boundaryFieldRef();
287
288 forAll(bpPrime, patchi)
289 {
290 if (!bpPrime[patchi].coupled())
291 {
292 bpPrime[patchi] == 0;
293 }
294 }
295
296 return tpPrime;
297}
298
299
300Foam::tmp<Foam::surfaceScalarField>
302{
303 return fvc::interpolate(pPrime());
304}
305
306
307Foam::tmp<Foam::volSymmTensorField>
309{
310 return devRhoReff(U_);
311}
312
313
314Foam::tmp<Foam::volSymmTensorField>
316(
317 const volVectorField& U
318) const
319{
321 (
322 IOobject::groupName("devRhoReff", U.group()),
324 (
325 - (rho_*nut_)
327 - ((rho_*lambda_)*fvc::div(phi_))*symmTensor::I
328 )
329 );
330}
331
332
333Foam::tmp<Foam::fvVectorMatrix>
335(
337) const
338{
339 return
340 (
341 - fvm::laplacian(rho_*nut_, U)
342 - fvc::div
343 (
344 (rho_*nut_)*dev2(T(fvc::grad(U)))
345 + ((rho_*lambda_)*fvc::div(phi_))
346 *dimensioned<symmTensor>("I", dimless, symmTensor::I)
347 )
348 );
349}
350
351
353{
354 // Local references
355 const twoPhaseSystem& fluid = refCast<const twoPhaseSystem>(phase_.fluid());
356 volScalarField alpha(max(alpha_, scalar(0)));
357 const volScalarField& rho = phase_.rho();
358 const surfaceScalarField& alphaRhoPhi = alphaRhoPhi_;
359 const volVectorField& U = U_;
360 const volVectorField& Uc_ = fluid.otherPhase(phase_).U();
361
362 const scalar sqrtPi = sqrt(constant::mathematical::pi);
363 dimensionedScalar ThetaSmall("ThetaSmall", Theta_.dimensions(), 1.0e-6);
364 dimensionedScalar ThetaSmallSqrt(sqrt(ThetaSmall));
365
366 tmp<volScalarField> tda(phase_.d());
367 const volScalarField& da = tda();
368
369 tmp<volTensorField> tgradU(fvc::grad(U_));
370 const volTensorField& gradU(tgradU());
371 volSymmTensorField D(symm(gradU));
372
373 // Calculating the radial distribution function
374 gs0_ = radialModel_->g0(alpha, alphaMinFriction_, alphaMax_);
375
376 if (!equilibrium_)
377 {
378 // Particle viscosity (Table 3.2, p.47)
379 nut_ = viscosityModel_->nu(alpha, Theta_, gs0_, rho, da, e_);
380
381 volScalarField ThetaSqrt("sqrtTheta", sqrt(Theta_));
382
383 // Bulk viscosity p. 45 (Lun et al. 1984).
384 lambda_ = (4.0/3.0)*sqr(alpha)*da*gs0_*(1.0 + e_)*ThetaSqrt/sqrtPi;
385
386 // Stress tensor, Definitions, Table 3.1, p. 43
388 (
389 rho*(2.0*nut_*D + (lambda_ - (2.0/3.0)*nut_)*tr(D)*I)
390 );
391
392 // Dissipation (Eq. 3.24, p.50)
393 volScalarField gammaCoeff
394 (
395 "gammaCoeff",
396 12.0*(1.0 - sqr(e_))
397 *max(sqr(alpha), residualAlpha_)
398 *rho*gs0_*(1.0/da)*ThetaSqrt/sqrtPi
399 );
400
401 // Drag
403 (
404 fluid.lookupSubModel<dragModel>
405 (
406 phase_,
407 fluid.otherPhase(phase_)
408 ).K()
409 );
410
411 // Eq. 3.25, p. 50 Js = J1 - J2
412 volScalarField J1("J1", 3.0*beta);
414 (
415 "J2",
416 0.25*sqr(beta)*da*magSqr(U - Uc_)
417 /(
418 max(alpha, residualAlpha_)*rho
419 *sqrtPi*(ThetaSqrt + ThetaSmallSqrt)
420 )
421 );
422
423 // particle pressure - coefficient in front of Theta (Eq. 3.22, p. 45)
424 volScalarField PsCoeff
425 (
426 granularPressureModel_->granularPressureCoeff
427 (
428 alpha,
429 gs0_,
430 rho,
431 e_
432 )
433 );
434
435 // 'thermal' conductivity (Table 3.3, p. 49)
436 kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rho, da, e_);
437
438 fv::options& fvOptions(fv::options::New(mesh_));
439
440 // Construct the granular temperature equation (Eq. 3.20, p. 44)
441 // NB. note that there are two typos in Eq. 3.20:
442 // Ps should be without grad
443 // the laplacian has the wrong sign
444 fvScalarMatrix ThetaEqn
445 (
446 1.5*
447 (
448 fvm::ddt(alpha, rho, Theta_)
449 + fvm::div(alphaRhoPhi, Theta_)
450 - fvc::Sp(fvc::ddt(alpha, rho) + fvc::div(alphaRhoPhi), Theta_)
451 )
452 - fvm::laplacian(kappa_, Theta_, "laplacian(kappa,Theta)")
453 ==
454 - fvm::SuSp((PsCoeff*I) && gradU, Theta_)
455 + (tau && gradU)
456 + fvm::Sp(-gammaCoeff, Theta_)
457 + fvm::Sp(-J1, Theta_)
458 + fvm::Sp(J2/(Theta_ + ThetaSmall), Theta_)
459 + fvOptions(alpha, rho, Theta_)
460 );
461
462 ThetaEqn.relax();
463 fvOptions.constrain(ThetaEqn);
464 ThetaEqn.solve();
465 fvOptions.correct(Theta_);
466 }
467 else
468 {
469 // Equilibrium => dissipation == production
470 // Eq. 4.14, p.82
471 volScalarField K1("K1", 2.0*(1.0 + e_)*rho*gs0_);
473 (
474 "K3",
475 0.5*da*rho*
476 (
477 (sqrtPi/(3.0*(3.0 - e_)))
478 *(1.0 + 0.4*(1.0 + e_)*(3.0*e_ - 1.0)*alpha*gs0_)
479 +1.6*alpha*gs0_*(1.0 + e_)/sqrtPi
480 )
481 );
482
484 (
485 "K2",
486 4.0*da*rho*(1.0 + e_)*alpha*gs0_/(3.0*sqrtPi) - 2.0*K3/3.0
487 );
488
489 volScalarField K4("K4", 12.0*(1.0 - sqr(e_))*rho*gs0_/(da*sqrtPi));
490
492 (
493 "trD",
494 alpha/(alpha + residualAlpha_)
495 *fvc::div(phi_)
496 );
497 volScalarField tr2D("tr2D", sqr(trD));
498 volScalarField trD2("trD2", tr(D & D));
499
500 volScalarField t1("t1", K1*alpha + rho);
501 volScalarField l1("l1", -t1*trD);
502 volScalarField l2("l2", sqr(t1)*tr2D);
504 (
505 "l3",
506 4.0
507 *K4
508 *alpha
509 *(2.0*K3*trD2 + K2*tr2D)
510 );
511
512 Theta_ = sqr
513 (
514 (l1 + sqrt(l2 + l3))
515 /(2.0*max(alpha, residualAlpha_)*K4)
516 );
517
518 kappa_ = conductivityModel_->kappa(alpha, Theta_, gs0_, rho, da, e_);
519 }
520
521 Theta_.clamp_range(0, 100);
522
523 {
524 // particle viscosity (Table 3.2, p.47)
525 nut_ = viscosityModel_->nu(alpha, Theta_, gs0_, rho, da, e_);
526
527 volScalarField ThetaSqrt("sqrtTheta", sqrt(Theta_));
528
529 // Bulk viscosity p. 45 (Lun et al. 1984).
530 lambda_ = (4.0/3.0)*sqr(alpha)*da*gs0_*(1.0 + e_)*ThetaSqrt/sqrtPi;
531
532 // Frictional pressure
534 (
535 frictionalStressModel_->frictionalPressure
536 (
537 phase_,
538 alphaMinFriction_,
539 alphaMax_
540 )
541 );
542
543 nuFric_ = frictionalStressModel_->nu
544 (
545 phase_,
546 alphaMinFriction_,
547 alphaMax_,
548 pf/rho,
549 D
550 );
551
552 // Limit viscosity and add frictional viscosity
553 nut_.clamp_max(maxNut_);
554
555 nuFric_ = min(nuFric_, maxNut_ - nut_);
556 nut_ += nuFric_;
557 }
558
559 if (debug)
560 {
561 Info<< typeName << ':' << nl
562 << " max(Theta) = " << max(Theta_).value() << nl
563 << " max(nut) = " << max(nut_).value() << endl;
564 }
565}
566
567
568// ************************************************************************* //
#define K3
Definition SHA1.C:145
#define K1
Definition SHA1.C:143
#define K4
Definition SHA1.C:146
#define K2
Definition SHA1.C:144
fv::options & fvOptions
twoPhaseSystem & fluid
static tmp< GeometricField< symmTensor, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< symmTensor >::calculatedType())
GeometricBoundaryField< scalar, fvPatchField, volMesh > Boundary
@ NO_REGISTER
Do not request registration (bool: false).
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
virtual tmp< volScalarField > epsilon() const
Return the turbulence kinetic energy dissipation rate.
virtual tmp< volSymmTensorField > R() const
Return the Reynolds stress tensor.
virtual void correct()
Solve the kinetic theory equations and correct the viscosity.
virtual tmp< volScalarField > k() const
Return the turbulence kinetic energy.
virtual tmp< fvVectorMatrix > divDevRhoReff(volVectorField &U) const
Return the source term for the momentum equation.
virtual tmp< volScalarField > pPrime() const
Return the phase-pressure'.
virtual tmp< surfaceScalarField > pPrimef() const
Return the face-phase-pressure'.
virtual tmp< volScalarField > omega() const
Return the specific dissipation rate.
virtual tmp< volSymmTensorField > devRhoReff() const
Return the effective stress tensor.
virtual bool read()
Re-read model coefficients if they have changed.
static const SymmTensor I
Definition SymmTensor.H:74
static options & New(const fvMesh &mesh)
Construct fvOptions and register to database if not present otherwise lookup and return.
Definition fvOptions.C:116
U
Definition pEqn.H:72
const volScalarField & T
bool coupled
dynamicFvMesh & mesh
#define NotImplemented
Issue a FatalErrorIn for a function not currently implemented.
Definition error.H:688
auto & name
word timeName
Definition getTimeIndex.H:3
List< ReturnType > get(const UPtrList< T > &list, const AccessOp &aop)
List of values generated by applying the access operation to each list item.
RASModel< EddyDiffusivity< turbulenceModel > > RASModel
constexpr scalar pi(M_PI)
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.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh > > grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcGrad.C:47
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition fvcDiv.C:42
tmp< GeometricField< Type, fvPatchField, volMesh > > ddt(const dimensioned< Type > dt, const fvMesh &mesh)
Definition fvcDdt.C:40
tmp< GeometricField< Type, fvPatchField, volMesh > > Sp(const volScalarField &sp, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvcSup.C:62
tmp< fvMatrix< Type > > laplacian(const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
tmp< fvMatrix< Type > > div(const surfaceScalarField &flux, const GeometricField< Type, fvPatchField, volMesh > &vf, const word &name)
Definition fvmDiv.C:41
tmp< fvMatrix< Type > > ddt(const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition fvmDdt.C:41
zeroField SuSp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
type
Types of root.
Definition Roots.H:53
dimensionedSymmTensor dev2(const dimensionedSymmTensor &dt)
const dimensionSet dimViscosity
Type & refCast(U &obj)
A dynamic_cast (for references) to Type reference.
Definition typeInfo.H:172
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:40
dimensionedSymmTensor symm(const dimensionedSymmTensor &dt)
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.
GeometricField< vector, fvPatchField, volMesh > volVectorField
const dimensionSet dimless
Dimensionless.
dimensionedScalar tr(const dimensionedSphericalTensor &dt)
bool read(const char *buf, int32_t &val)
Same as readInt32.
Definition int32.H:127
dimensionedSymmTensor sqr(const dimensionedVector &dv)
fvMatrix< scalar > fvScalarMatrix
GeometricField< scalar, fvPatchField, volMesh > volScalarField
SymmTensor< Cmpt > devTwoSymm(const SymmTensor< Cmpt > &st)
Return the deviatoric part of twice the symmetric part of a SymmTensor.
messageStream Info
Information stream (stdout output on master, null elsewhere).
static const Identity< scalar > I
Definition Identity.H:100
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition Ostream.H:519
GeometricField< tensor, fvPatchField, volMesh > volTensorField
dimensionedScalar sqrt(const dimensionedScalar &ds)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition hashSets.C:26
ThermalDiffusivity< PhaseCompressibleTurbulenceModel< phaseModel > > phaseCompressibleTurbulenceModel
Typedef for phaseCompressibleTurbulenceModel.
GeometricField< symmTensor, fvPatchField, volMesh > volSymmTensorField
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dimensioned< typename typeOfMag< Type >::type > magSqr(const dimensioned< Type > &dt)
constexpr char nl
The newline '\n' character (0x0a).
Definition Ostream.H:50
volScalarField & alpha
const dimensionedScalar & D
dimensionedScalar beta("beta", dimless/dimTemperature, laminarTransport)
#define forAll(list, i)
Loop across all elements in list.
Definition stdFoam.H:299