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-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 
26 #include "LehrMilliesMewes.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace binaryBreakupModels
37 {
40  (
44  );
45 }
46 }
47 }
48 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
54 (
55  const populationBalanceModel& popBal,
56  const dictionary& dict
57 )
58 :
59  binaryBreakupModel(popBal, dict)
60 {}
61 
62 
63 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
64 
67 (
68  volScalarField::Internal& binaryBreakupRate,
69  const label i,
70  const label j
71 )
72 {
73  const sizeGroup& fi = popBal_.sizeGroups()[i];
74  const sizeGroup& fj = popBal_.sizeGroups()[j];
75 
76  const volScalarField::Internal& rhoc = popBal_.continuousPhase().rho();
77 
78  tmp<volScalarField> tsigma(popBal_.sigmaWithContinuousPhase(fj.phase()));
79  const volScalarField::Internal& sigma = tsigma();
80 
81  tmp<volScalarField> tepsilonc(popBal_.continuousTurbulence().epsilon());
82  const volScalarField::Internal& epsilonc = tepsilonc();
83 
85  (
86  pow(sigma/rhoc, 3.0/5.0)/pow(epsilonc, 2.0/5.0)
87  );
88 
89  // Reset of dimension to pure length to avoid problems in transcendental
90  // functions due to small exponents
92 
94  (
95  pow(sigma/rhoc, 2.0/5.0)/pow(epsilonc, 3.0/5.0)
96  );
97 
98  binaryBreakupRate +=
99  0.5*pow(fj.dSph()/L, 5.0/3.0)
100  *exp(-sqrt(2.0)/pow3(fj.dSph()/L))
101  *6/pow(pi, 1.5)/pow3(fi.dSph()/L)
102  *exp(-9.0/4.0*sqr(log(pow(2.0, 0.4)*fi.dSph()/L)))
103  /max(1 + erf(1.5*log(pow(2.0, 1.0/15.0)*fj.dSph()/L)), small)
104  /(T*pow3(L));
105 }
106 
107 
108 // ************************************************************************* //
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...
const dimensionSet & dimensions() const
Return dimensions.
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::Internal &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: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
void reset(const dimensionSet &)
Definition: dimensionSet.C:126
A class for managing temporary objects.
Definition: tmp.H:55
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m^2/K^4].
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
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
dimensionedScalar erf(const dimensionedScalar &ds)
const dimensionSet dimLength
dimensionedScalar log(const dimensionedScalar &ds)
void pow(LagrangianPatchField< typename powProduct< Type, r >::type > &f, const LagrangianPatchField< Type > &f1)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void pow3(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict