LehrMilliesMewesCoalescence.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 
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace diameterModels
34 {
35 namespace coalescenceModels
36 {
39  (
43  );
44 }
45 }
46 }
47 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
54 (
55  const populationBalanceModel& popBal,
56  const dictionary& dict
57 )
58 :
59  coalescenceModel(popBal, dict),
60  uCrit_
61  (
62  dimensionedScalar::lookupOrDefault("uCrit", dict, dimVelocity, 0.08)
63  ),
64  alphaMax_
65  (
66  dimensionedScalar::lookupOrDefault("alphaMax", dict, dimless, 0.6)
67  )
68 {}
69 
70 
71 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
72 
73 void
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 
85  const volScalarField uChar
86  (
87  max
88  (
89  sqrt(2.0)*cbrt(popBal_.continuousTurbulence().epsilon())
90  *sqrt(cbrt(sqr(fi.dSph())) + cbrt(sqr(fj.dSph()))),
91  mag(fi.phase().U() - fj.phase().U())
92  )
93  );
94 
95  coalescenceRate +=
96  pi/4*sqr(fi.dSph() + fj.dSph())*min(uChar, uCrit_)
97  *exp
98  (
99  - sqr(cbrt(alphaMax_)
100  /cbrt(max(popBal_.alphas(), fi.phase().residualAlpha())) - 1)
101  );
102 }
103 
104 
105 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Base class for coalescence models.
Model of Lehr et al. (2002). The coalescence rate is calculated by.
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
LehrMilliesMewesCoalescence(const populationBalanceModel &popBal, const dictionary &dict)
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:50
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:36
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:160
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:133
virtual tmp< volVectorField > U() const =0
Return the velocity.
addToRunTimeSelectionTable(coalescenceModel, AdachiStuartFokkink, dictionary)
defineTypeNameAndDebug(AdachiStuartFokkink, 0)
Namespace for OpenFOAM.
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)
const dimensionSet dimless
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
const dimensionSet dimVelocity
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict