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 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 #include "phaseCompressibleTurbulenceModel.H"
29 #include "mathematicalConstants.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace binaryBreakupModels
38 {
39  defineTypeNameAndDebug(LehrMilliesMewes, 0);
41  (
42  binaryBreakupModel,
43  LehrMilliesMewes,
44  dictionary
45  );
46 }
47 }
48 }
49 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
55 (
56  const populationBalanceModel& popBal,
57  const dictionary& dict
58 )
59 :
60  binaryBreakupModel(popBal, dict)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
68 (
69  volScalarField& binaryBreakupRate,
70  const label i,
71  const label j
72 )
73 {
74  const phaseModel& continuousPhase = popBal_.continuousPhase();
75  const sizeGroup& fi = popBal_.sizeGroups()[i];
76  const sizeGroup& fj = popBal_.sizeGroups()[j];
77 
79  (
80  pow
81  (
82  popBal_.sigmaWithContinuousPhase(fj.phase())/continuousPhase.rho(),
83  3.0/5.0
84  )
86  );
87 
88  // Reset of dimension to pure length to avoid problems in transcendental
89  // functions due to small exponents
91 
92  const volScalarField T
93  (
94  pow
95  (
96  popBal_.sigmaWithContinuousPhase(fj.phase())/continuousPhase.rho(),
97  2.0/5.0
98  )
100  );
101 
102  binaryBreakupRate +=
103  0.5*pow(fj.d()/L, 5.0/3.0)
104  *exp(-sqrt(2.0)/pow3(fj.d()/L))
105  *6.0/pow(pi, 1.5)/pow3(fi.d()/L)
106  *exp(-9.0/4.0*sqr(log(pow(2.0, 0.4)*fi.d()/L)))
107  /max(1.0 + erf(1.5*log(pow(2.0, 1.0/15.0)*fj.d()/L)), SMALL)
108  /(T*pow3(L));
109 }
110 
111 
112 // ************************************************************************* //
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
LehrMilliesMewes(const populationBalanceModel &popBal, const dictionary &dict)
dimensionedScalar log(const dimensionedScalar &ds)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar sqrt(const dimensionedScalar &ds)
static const scalar SMALL
Definition: scalar.H:112
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const dimensionSet & dimensions() const
Return dimensions.
dimensionedScalar exp(const dimensionedScalar &ds)
const phaseModel & continuousPhase() const
Return continuous phase.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
dimensionedScalar erf(const dimensionedScalar &ds)
defineTypeNameAndDebug(combustionModel, 0)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
void reset(const dimensionSet &)
Definition: dimensionSet.C:108
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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.