LehrMilliesMewes.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2019 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "LehrMilliesMewes.H"
30#include "phaseCompressibleTurbulenceModel.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34
35namespace Foam
36{
37namespace diameterModels
38{
39namespace binaryBreakupModels
40{
43 (
47 );
48}
49}
50}
51
53
54// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
55
57(
58 const populationBalanceModel& popBal,
59 const dictionary& dict
60)
61:
63{}
64
65
66// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
67
70(
71 volScalarField& binaryBreakupRate,
72 const label i,
73 const label j
74)
75{
76 const phaseModel& continuousPhase = popBal_.continuousPhase();
77 const sizeGroup& fi = popBal_.sizeGroups()[i];
78 const sizeGroup& fj = popBal_.sizeGroups()[j];
79
81 (
82 pow
83 (
84 popBal_.sigmaWithContinuousPhase(fj.phase())/continuousPhase.rho(),
85 3.0/5.0
86 )
87 /pow(popBal_.continuousTurbulence().epsilon(), 2.0/5.0)
88 );
89
90 // Reset of dimension to pure length to avoid problems in transcendental
91 // functions due to small exponents
92 L.dimensions().reset(dimLength);
93
94 const volScalarField T
95 (
96 pow
97 (
98 popBal_.sigmaWithContinuousPhase(fj.phase())/continuousPhase.rho(),
99 2.0/5.0
100 )
101 /pow(popBal_.continuousTurbulence().epsilon(), 3.0/5.0)
102 );
103
104 binaryBreakupRate +=
105 0.5*pow(fj.d()/L, 5.0/3.0)
106 *exp(-sqrt(2.0)/pow3(fj.d()/L))
107 *6.0/pow(pi, 1.5)/pow3(fi.d()/L)
108 *exp(-9.0/4.0*sqr(log(pow(2.0, 0.4)*fi.d()/L)))
109 /max(1.0 + erf(1.5*log(pow(2.0, 1.0/15.0)*fj.d()/L)), SMALL)
110 /(T*pow3(L));
111}
112
113
114// ************************************************************************* //
Macros for easy insertion into run-time selection tables.
#define addToRunTimeSelectionTable(baseType, thisType, argNames)
Add to construction table with typeName as the key.
Base class for binary breakup models which give the breakup rate between a sizeGroup pair directly,...
Model of Lehr et al. (2002). The breakup rate is calculated by.
virtual void addToBinaryBreakupRate(volScalarField &binaryBreakupRate, const label i, const label j)
Add to binary breakupRate.
Class that solves the univariate population balance equation by means of a class method (also called ...
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:99
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroupI.H:52
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:38
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:61
const dimensionedScalar & rho() const
Definition: phaseModel.H:171
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
const volScalarField & T
constexpr scalar pi(M_PI)
Namespace for OpenFOAM.
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:52
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar erf(const dimensionedScalar &ds)
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dictionary dict
const vector L(dict.get< vector >("L"))