Luo.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) 2019 OpenFOAM Foundation
9  Copyright (C) 2020 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 "Luo.H"
31 #include "mathematicalConstants.H"
32 #include "phaseCompressibleTurbulenceModel.H"
33 #include "virtualMassModel.H"
34 #include "phaseSystem.H"
35 
36 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
37 
38 namespace Foam
39 {
40 namespace diameterModels
41 {
42 namespace coalescenceModels
43 {
44  defineTypeNameAndDebug(Luo, 0);
46  (
47  coalescenceModel,
48  Luo,
49  dictionary
50  );
51 }
52 }
53 }
54 
56 
57 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58 
60 Luo
61 (
62  const populationBalanceModel& popBal,
63  const dictionary& dict
64 )
65 :
66  coalescenceModel(popBal, dict),
67  beta_(dimensionedScalar::getOrDefault("beta", dict, dimless, 2.05)),
68  C1_(dimensionedScalar::getOrDefault("C1", dict, dimless, 1.0))
69 {}
70 
71 
72 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73 
76 (
77  volScalarField& coalescenceRate,
78  const label i,
79  const label j
80 )
81 {
82  const sizeGroup& fi = popBal_.sizeGroups()[i];
83  const sizeGroup& fj = popBal_.sizeGroups()[j];
84  const phaseModel& continuousPhase = popBal_.continuousPhase();
85 
86  if
87  (
88  popBal_.fluid().foundSubModel<virtualMassModel>
89  (
90  fi.phase(),
91  popBal_.continuousPhase()
92  )
93  )
94  {
95  const virtualMassModel& vm =
96  popBal_.fluid().lookupSubModel<virtualMassModel>
97  (
98  fi.phase(),
99  popBal_.continuousPhase()
100  );
101 
102  const dimensionedScalar xi = fi.d()/fj.d();
103 
104  const volScalarField uij
105  (
106  sqrt(beta_)
107  *cbrt(popBal_.continuousTurbulence().epsilon()*fi.d())
108  *sqrt(1.0 + pow(xi, -2.0/3.0))
109  );
110 
111  coalescenceRate +=
112  pi/4.0*sqr(fi.d() + fj.d())*uij
113  *exp
114  (
115  - C1_
116  *sqrt(0.75*(1.0 + sqr(xi))*(1.0 + pow3(xi)))
117  /(
118  sqrt(fi.phase().rho()/continuousPhase.rho()
119  + vm.Cvm())*pow3(1.0 + xi)
120  )
121  *sqrt
122  (
123  continuousPhase.rho()*fi.d()*sqr(uij)
124  /popBal_.sigmaWithContinuousPhase(fi.phase())
125  )
126  );
127  }
128  else
129  {
131  << "A virtual mass model for " << fi.phase().name() << " in "
132  << popBal_.continuousPhase().name() << " is not specified. This is "
133  << "required by the Luo coalescence model." << exit(FatalError);
134  }
135 }
136 
137 
138 // ************************************************************************* //
Foam::virtualMassModel
Definition: virtualMassModel.H:55
Foam::diameterModels::coalescenceModel
Base class for coalescence models.
Definition: coalescenceModel.H:52
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
mathematicalConstants.H
Foam::diameterModels::sizeGroup::phase
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:38
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::diameterModels::coalescenceModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(coalescenceModel, constantCoalescence, dictionary)
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::virtualMassModel::Cvm
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::diameterModels::sizeGroup
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:96
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::diameterModels::populationBalanceModel
Class that solves the univariate population balance equation by means of a class method (also called ...
Definition: populationBalanceModel.H:179
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam::dimensioned< scalar >
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
Foam::phaseModel::name
const word & name() const
Definition: phaseModel.H:140
Foam::diameterModels::sizeGroup::d
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroupI.H:52
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Luo.H
Foam::diameterModels::coalescenceModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constantCoalescence, 0)
Foam::phaseModel::rho
const dimensionedScalar & rho() const
Definition: phaseModel.H:167
Foam::diameterModels::coalescenceModels::Luo::addToCoalescenceRate
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition: Luo.C:76
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::diameterModels::coalescenceModels::Luo::Luo
Luo(const populationBalanceModel &popBal, const dictionary &dict)
Definition: Luo.C:61
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189