IshiiZuber.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) 2014-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 "IshiiZuber.H"
27 #include "phasePair.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace dragModels
35 {
36  defineTypeNameAndDebug(IshiiZuber, 0);
37  addToRunTimeSelectionTable(dragModel, IshiiZuber, dictionary);
38 }
39 }
40 
41 
42 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
43 
45 (
46  const dictionary& dict,
47  const phasePair& pair,
48  const bool registerObject
49 )
50 :
51  dragModel(dict, pair, registerObject)
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
56 
58 {}
59 
60 
61 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
62 
65 {
66  const volScalarField Re(pair_.Re());
67  const volScalarField Eo(pair_.Eo());
68 
69  const volScalarField mud(pair_.dispersed().thermo().mu());
70  const volScalarField muc(pair_.continuous().thermo().mu());
71 
72  const volScalarField muStar((mud + 0.4*muc)/(mud + muc));
73 
74  const volScalarField muMix
75  (
76  muc*pow(max(1 - pair_.dispersed(), scalar(1e-3)), -2.5*muStar)
77  );
78 
79  const volScalarField ReM(Re*muc/muMix);
80  const volScalarField CdRe
81  (
82  pos0(1000 - ReM)*24*(1 + 0.1*pow(ReM, 0.75))
83  + neg(1000 - ReM)*0.44*ReM
84  );
85 
86  volScalarField F((muc/muMix)*sqrt(1 - pair_.dispersed()));
87  F.max(1e-3);
88 
89  const volScalarField Ealpha((1 + 17.67*pow(F, 0.8571428))/(18.67*F));
90 
91  const volScalarField CdReEllipse(Ealpha*0.6666*sqrt(Eo)*Re);
92 
93  return
94  pos0(CdReEllipse - CdRe)
95  *min(CdReEllipse, Re*sqr(1 - pair_.dispersed())*2.66667)
96  + neg(CdReEllipse - CdRe)*CdRe;
97 }
98 
99 
100 // ************************************************************************* //
const dimensionedScalar & F
Faraday constant: default SI units: [C/mol].
dictionary dict
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
virtual ~IshiiZuber()
Destructor.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
tmp< volScalarField > Eo() const
Eotvos number.
virtual const phaseModel & continuous() const
Continuous phase.
virtual tmp< volScalarField > mu() const
Dynamic viscosity of mixture [kg/m/s].
Definition: rhoThermo.C:188
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:120
tmp< volScalarField > Re() const
Reynolds number.
IshiiZuber(const dictionary &dict, const phasePair &pair, const bool registerObject)
Construct from a dictionary and a phase pair.
virtual tmp< volScalarField > CdRe() const
Drag coefficient.
virtual const phaseModel & dispersed() const
Dispersed phase.
const phasePair & pair_
Phase pair.
Definition: dragModel.H:62
dimensionedScalar pos0(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const doubleScalar e
Elementary charge.
Definition: doubleScalar.H:105
A class for managing temporary objects.
Definition: PtrList.H:53
scalarField Re(const UList< complex > &cf)
Definition: complexFields.C:97
Namespace for OpenFOAM.