LehrMilliesMewesCoalescence.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 
28 #include "mathematicalConstants.H"
29 #include "phaseCompressibleTurbulenceModel.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace coalescenceModels
38 {
39  defineTypeNameAndDebug(LehrMilliesMewesCoalescence, 0);
41  (
42  coalescenceModel,
43  LehrMilliesMewesCoalescence,
44  dictionary
45  );
46 }
47 }
48 }
49 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  coalescenceModel(popBal, dict),
62  uCrit_
63  (
64  dimensionedScalar::lookupOrDefault("uCrit", dict, dimVelocity, 0.08)
65  ),
66  alphaMax_
67  (
68  dimensionedScalar::lookupOrDefault("alphaMax", dict, dimless, 0.6)
69  )
70 {}
71 
72 
73 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
74 
75 void
78 (
79  volScalarField& coalescenceRate,
80  const label i,
81  const label j
82 )
83 {
84  const sizeGroup& fi = popBal_.sizeGroups()[i];
85  const sizeGroup& fj = popBal_.sizeGroups()[j];
86 
87  const volScalarField uChar
88  (
89  max
90  (
92  *sqrt(cbrt(sqr(fi.d())) + cbrt(sqr(fj.d()))),
93  mag(fi.phase().U() - fj.phase().U())
94  )
95  );
96 
97  coalescenceRate +=
98  pi/4.0*sqr(fi.d() + fj.d())*min(uChar, uCrit_)
99  *exp
100  (
101  - sqr(cbrt(alphaMax_)
102  /cbrt(max(popBal_.alphas(), fi.phase().residualAlpha())) - 1.0)
103  );
104 }
105 
106 
107 // ************************************************************************* //
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
dimensionedScalar sqrt(const dimensionedScalar &ds)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
dimensionedScalar exp(const dimensionedScalar &ds)
LehrMilliesMewesCoalescence(const populationBalanceModel &popBal, const dictionary &dict)
dimensionedScalar cbrt(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
const volScalarField & alphas() const
Return total void of phases belonging to this populationBalance.
dimensioned< scalar > mag(const dimensioned< Type > &)
Namespace for OpenFOAM.
const dimensionSet dimVelocity