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-2025 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 
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace coalescenceModels
37 {
40  (
44  );
45 }
46 }
47 }
48 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const populationBalanceModel& popBal,
57  const dictionary& dict
58 )
59 :
60  coalescenceModel(popBal, dict),
61  uCrit_("uCrit", dimVelocity, dict, 0.08),
62  alphaMax_("alphaMax", dimless, dict, 0.6)
63 {}
64 
65 
66 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67 
70 (
71  volScalarField::Internal& coalescenceRate,
72  const label i,
73  const label j
74 )
75 {
76  const sizeGroup& fi = popBal_.sizeGroups()[i];
77  const sizeGroup& fj = popBal_.sizeGroups()[j];
78 
79  tmp<volScalarField> tepsilonc(popBal_.continuousTurbulence().epsilon());
80  const volScalarField::Internal& epsilonc = tepsilonc();
81 
82  const volScalarField::Internal uChar
83  (
84  max
85  (
86  sqrt(2.0)
87  *cbrt(epsilonc)
88  *sqrt(cbrt(sqr(fi.dSph())) + cbrt(sqr(fj.dSph()))),
89  mag(fi.phase().U()()() - fj.phase().U()()())
90  )
91  );
92 
93  coalescenceRate +=
94  pi/4
95  *sqr(fi.dSph() + fj.dSph())
96  *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.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Base class for coalescence models.
Model of Lehr et al. (2002). The coalescence rate is calculated by.
virtual void addToCoalescenceRate(volScalarField::Internal &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:96
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 keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:169
virtual tmp< volVectorField > U() const =0
Return the velocity.
A class for managing temporary objects.
Definition: tmp.H:55
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
const dimensionSet dimless
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void cbrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimVelocity
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict