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 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"
28 #include "mathematicalConstants.H"
29 #include "phaseCompressibleTurbulenceModel.H"
30 #include "virtualMassModel.H"
31 #include "phaseSystem.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace diameterModels
38 {
39 namespace coalescenceModels
40 {
41  defineTypeNameAndDebug(Luo, 0);
43  (
44  coalescenceModel,
45  Luo,
46  dictionary
47  );
48 }
49 }
50 }
51 
53 
54 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55 
57 Luo
58 (
59  const populationBalanceModel& popBal,
60  const dictionary& dict
61 )
62 :
63  coalescenceModel(popBal, dict),
64  beta_(dimensionedScalar::lookupOrDefault("beta", dict, dimless, 2.05)),
65  C1_(dimensionedScalar::lookupOrDefault("C1", dict, dimless, 1.0))
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
73 (
74  volScalarField& coalescenceRate,
75  const label i,
76  const label j
77 )
78 {
79  const sizeGroup& fi = popBal_.sizeGroups()[i];
80  const sizeGroup& fj = popBal_.sizeGroups()[j];
81  const phaseModel& continuousPhase = popBal_.continuousPhase();
82 
83  if
84  (
85  popBal_.fluid().foundSubModel<virtualMassModel>
86  (
87  fi.phase(),
89  )
90  )
91  {
92  const virtualMassModel& vm =
93  popBal_.fluid().lookupSubModel<virtualMassModel>
94  (
95  fi.phase(),
97  );
98 
99  const dimensionedScalar xi = fi.d()/fj.d();
100 
101  const volScalarField uij
102  (
103  sqrt(beta_)
105  *sqrt(1.0 + pow(xi, -2.0/3.0))
106  );
107 
108  coalescenceRate +=
109  pi/4.0*sqr(fi.d() + fj.d())*uij
110  *exp
111  (
112  - C1_
113  *sqrt(0.75*(1.0 + sqr(xi))*(1.0 + pow3(xi)))
114  /(
115  sqrt(fi.phase().rho()/continuousPhase.rho()
116  + vm.Cvm())*pow3(1.0 + xi)
117  )
118  *sqrt
119  (
120  continuousPhase.rho()*fi.d()*sqr(uij)
121  /popBal_.sigmaWithContinuousPhase(fi.phase())
122  )
123  );
124  }
125  else
126  {
128  << "A virtual mass model for " << fi.phase().name() << " in "
129  << popBal_.continuousPhase().name() << " is not specified. This is "
130  << "required by the Luo coalescence model." << exit(FatalError);
131  }
132 }
133 
134 
135 // ************************************************************************* //
Luo(const populationBalanceModel &popBal, const dictionary &dict)
dictionary dict
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
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:319
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const word & name() const
Definition: phaseModel.H:109
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:52
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 foundSubModel(const phasePair &key) const
Check availability of a sub model for a given phase pair.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
Namespace for OpenFOAM.