LehrMilliesMewes.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 "LehrMilliesMewes.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace diameterModels
34 {
35 namespace binaryBreakupModels
36 {
39  (
43  );
44 }
45 }
46 }
47 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
53 (
54  const populationBalanceModel& popBal,
55  const dictionary& dict
56 )
57 :
58  binaryBreakupModel(popBal, dict)
59 {}
60 
61 
62 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
63 
66 (
67  volScalarField& binaryBreakupRate,
68  const label i,
69  const label j
70 )
71 {
72  const phaseModel& continuousPhase = popBal_.continuousPhase();
73  const sizeGroup& fi = popBal_.sizeGroups()[i];
74  const sizeGroup& fj = popBal_.sizeGroups()[j];
75 
77  (
78  pow
79  (
80  popBal_.sigmaWithContinuousPhase(fj.phase())/continuousPhase.rho(),
81  3.0/5.0
82  )
83  /pow(popBal_.continuousTurbulence().epsilon(), 2.0/5.0)
84  );
85 
86  // Reset of dimension to pure length to avoid problems in transcendental
87  // functions due to small exponents
89 
90  const volScalarField T
91  (
92  pow
93  (
94  popBal_.sigmaWithContinuousPhase(fj.phase())/continuousPhase.rho(),
95  2.0/5.0
96  )
97  /pow(popBal_.continuousTurbulence().epsilon(), 3.0/5.0)
98  );
99 
100  binaryBreakupRate +=
101  0.5*pow(fj.dSph()/L, 5.0/3.0)
102  *exp(-sqrt(2.0)/pow3(fj.dSph()/L))
103  *6/pow(pi, 1.5)/pow3(fi.dSph()/L)
104  *exp(-9.0/4.0*sqr(log(pow(2.0, 0.4)*fi.dSph()/L)))
105  /max(1 + erf(1.5*log(pow(2.0, 1.0/15.0)*fj.dSph()/L)), small)
106  /(T*pow3(L));
107 }
108 
109 
110 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
const dimensionSet & dimensions() const
Return dimensions.
Generic GeometricField class.
Base class for binary breakup models that provide a breakup rate between a size class pair directly,...
Model of Lehr et al. (2002). The breakup rate is calculated by.
LehrMilliesMewes(const populationBalanceModel &popBal, const dictionary &dict)
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
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
void reset(const dimensionSet &)
Definition: dimensionSet.C:108
virtual const volScalarField & rho() const =0
Return the density field.
addToRunTimeSelectionTable(binaryBreakupModel, LehrMilliesMewes, dictionary)
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)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
const dimensionSet dimLength
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict