All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Luo.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) 2019-2022 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 "Luo.H"
29 #include "phaseSystem.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace coalescenceModels
38 {
39  defineTypeNameAndDebug(Luo, 0);
41  (
42  coalescenceModel,
43  Luo,
44  dictionary
45  );
46 }
47 }
48 }
49 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 Luo
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  coalescenceModel(popBal, dict),
62  beta_(dimensionedScalar::lookupOrDefault("beta", dict, dimless, 2.05)),
63  C1_(dimensionedScalar::lookupOrDefault("C1", dict, dimless, 1))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
71 (
72  volScalarField& coalescenceRate,
73  const label i,
74  const label j
75 )
76 {
77  const sizeGroup& fi = popBal_.sizeGroups()[i];
78  const sizeGroup& fj = popBal_.sizeGroups()[j];
79  const phaseModel& continuousPhase = popBal_.continuousPhase();
80 
81  if
82  (
84  <virtualMassModels::dispersedVirtualMassModel>
85  (
86  dispersedPhaseInterface(fi.phase(), popBal_.continuousPhase())
87  )
88  )
89  {
90  const virtualMassModels::dispersedVirtualMassModel& vm =
92  <virtualMassModels::dispersedVirtualMassModel>
93  (
94  dispersedPhaseInterface(fi.phase(), popBal_.continuousPhase())
95  );
96 
97  const dimensionedScalar xi = fi.dSph()/fj.dSph();
98 
99  const volScalarField uij
100  (
101  sqrt(beta_)
102  *cbrt(popBal_.continuousTurbulence().epsilon()*fi.dSph())
103  *sqrt(1 + pow(xi, -2.0/3.0))
104  );
105 
106  coalescenceRate +=
107  pi/4*sqr(fi.dSph() + fj.dSph())*uij
108  *exp
109  (
110  - C1_
111  *sqrt(0.75*(1 + sqr(xi))*(1 + pow3(xi)))
112  /(
113  sqrt(fi.phase().rho()/continuousPhase.rho()
114  + vm.Cvm())*pow3(1 + xi)
115  )
116  *sqrt
117  (
118  continuousPhase.rho()*fi.dSph()*sqr(uij)
119  /popBal_.sigmaWithContinuousPhase(fi.phase())
120  )
121  );
122  }
123  else
124  {
126  << "A virtual mass model for " << fi.phase().name() << " in "
127  << popBal_.continuousPhase().name() << " is not specified. This is "
128  << "required by the Luo coalescence model." << exit(FatalError);
129  }
130 }
131 
132 
133 // ************************************************************************* //
Luo(const populationBalanceModel &popBal, const dictionary &dict)
dictionary dict
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:306
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const phaseCompressible::momentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const word & name() const
Definition: phaseModel.H:110
const dimensionSet dimless
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
dimensionedScalar exp(const dimensionedScalar &ds)
const phaseModel & continuousPhase() const
Return continuous phase.
const phaseSystem & fluid() const
Return reference to the phaseSystem.
dimensionedScalar cbrt(const dimensionedScalar &ds)
bool foundInterfacialModel(const phaseInterface &interface) const
Check availability of a sub model for a given interface.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const ModelType & lookupInterfacialModel(const phaseInterface &interface) const
Return a sub model for an interface.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Namespace for OpenFOAM.