LehrMilliesMewesCoalescence.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 Copyright (C) 2020 OpenCFD Ltd.
10-------------------------------------------------------------------------------
11License
12 This file is part of OpenFOAM.
13
14 OpenFOAM is free software: you can redistribute it and/or modify it
15 under the terms of the GNU General Public License as published by
16 the Free Software Foundation, either version 3 of the License, or
17 (at your option) any later version.
18
19 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22 for more details.
23
24 You should have received a copy of the GNU General Public License
25 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26
27\*---------------------------------------------------------------------------*/
28
32#include "phaseCompressibleTurbulenceModel.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35
36namespace Foam
37{
38namespace diameterModels
39{
40namespace coalescenceModels
41{
44 (
48 );
49}
50}
51}
52
54
55// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56
57Foam::diameterModels::coalescenceModels::LehrMilliesMewesCoalescence::
58LehrMilliesMewesCoalescence
59(
60 const populationBalanceModel& popBal,
61 const dictionary& dict
62)
63:
64 coalescenceModel(popBal, dict),
65 uCrit_
66 (
67 dimensionedScalar::getOrDefault("uCrit", dict, dimVelocity, 0.08)
68 ),
69 alphaMax_
70 (
71 dimensionedScalar::getOrDefault("alphaMax", dict, dimless, 0.6)
72 )
73{}
74
75
76// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
77
78void
81(
82 volScalarField& coalescenceRate,
83 const label i,
84 const label j
85)
86{
87 const sizeGroup& fi = popBal_.sizeGroups()[i];
88 const sizeGroup& fj = popBal_.sizeGroups()[j];
89
90 const volScalarField uChar
91 (
92 max
93 (
94 sqrt(2.0)*cbrt(popBal_.continuousTurbulence().epsilon())
95 *sqrt(cbrt(sqr(fi.d())) + cbrt(sqr(fj.d()))),
96 mag(fi.phase().U() - fj.phase().U())
97 )
98 );
99
100 coalescenceRate +=
101 pi/4.0*sqr(fi.d() + fj.d())*min(uChar, uCrit_)
102 *exp
103 (
104 - sqr(cbrt(alphaMax_)
105 /cbrt(max(popBal_.alphas(), fi.phase().residualAlpha())) - 1.0)
106 );
107}
108
109
110// ************************************************************************* //
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 coalescence models.
Model of Lehr et al. (2002). The coalescence rate is calculated by.
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
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
const volVectorField & U() const
Definition: phaseModel.H:176
const dimensionedScalar & residualAlpha() const
Return the residual phase-fraction for given phase.
Definition: phaseModel.C:154
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition: className.H:121
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)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimVelocity
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
dimensionedScalar cbrt(const dimensionedScalar &ds)
dictionary dict