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-2023 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"
27 #include "phaseSystem.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38 namespace coalescenceModels
39 {
42  (
44  Luo,
46  );
47 }
48 }
49 }
50 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
56 Luo
57 (
58  const populationBalanceModel& popBal,
59  const dictionary& dict
60 )
61 :
62  coalescenceModel(popBal, dict),
63  beta_(dimensionedScalar::lookupOrDefault("beta", dict, dimless, 2.05)),
64  C1_(dimensionedScalar::lookupOrDefault("C1", dict, dimless, 1))
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
72 (
73  volScalarField& coalescenceRate,
74  const label i,
75  const label j
76 )
77 {
78  const sizeGroup& fi = popBal_.sizeGroups()[i];
79  const sizeGroup& fj = popBal_.sizeGroups()[j];
80  const phaseModel& continuousPhase = popBal_.continuousPhase();
81 
82  if
83  (
84  popBal_.fluid().foundInterfacialModel
86  (
87  dispersedPhaseInterface(fi.phase(), popBal_.continuousPhase())
88  )
89  )
90  {
92  popBal_.fluid().lookupInterfacialModel
94  (
95  dispersedPhaseInterface(fi.phase(), popBal_.continuousPhase())
96  );
97 
98  const dimensionedScalar xi = fi.dSph()/fj.dSph();
99 
100  const volScalarField uij
101  (
102  sqrt(beta_)
103  *cbrt(popBal_.continuousTurbulence().epsilon()*fi.dSph())
104  *sqrt(1 + pow(xi, -2.0/3.0))
105  );
106 
107  coalescenceRate +=
108  pi/4*sqr(fi.dSph() + fj.dSph())*uij
109  *exp
110  (
111  - C1_
112  *sqrt(0.75*(1 + sqr(xi))*(1 + pow3(xi)))
113  /(
114  sqrt(fi.phase().rho()/continuousPhase.rho()
115  + vm.Cvm())*pow3(1 + xi)
116  )
117  *sqrt
118  (
119  continuousPhase.rho()*fi.dSph()*sqr(uij)
120  /popBal_.sigmaWithContinuousPhase(fi.phase())
121  )
122  );
123  }
124  else
125  {
127  << "A virtual mass model for " << fi.phase().name() << " in "
128  << popBal_.continuousPhase().name() << " is not specified. This is "
129  << "required by the Luo coalescence model." << exit(FatalError);
130  }
131 }
132 
133 
134 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Base class for coalescence models.
Model of Luo (1993). The coalescence rate is calculated by.
Definition: Luo.H:169
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition: Luo.C:72
Luo(const populationBalanceModel &popBal, const dictionary &dict)
Definition: Luo.C:57
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & dSph() const
Return representative spherical diameter of the sizeGroup.
Definition: sizeGroupI.H:51
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
Class to represent a interface between phases where one phase is considered dispersed within the othe...
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
virtual const volScalarField & rho() const =0
Return the density field.
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
addToRunTimeSelectionTable(coalescenceModel, AdachiStuartFokkink, dictionary)
defineTypeNameAndDebug(AdachiStuartFokkink, 0)
Namespace for OpenFOAM.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
error FatalError
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict