All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
CoulaloglouTavlaridesCoalescence.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) 2018-2021 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 
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace coalescenceModels
37 {
38  defineTypeNameAndDebug(CoulaloglouTavlaridesCoalescence, 0);
40  (
41  coalescenceModel,
42  CoulaloglouTavlaridesCoalescence,
43  dictionary
44  );
45 }
46 }
47 }
48 
49 
50 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
51 
54 (
55  const populationBalanceModel& popBal,
56  const dictionary& dict
57 )
58 :
59  coalescenceModel(popBal, dict),
60  C1_(dimensionedScalar::lookupOrDefault("C1", dict, dimless, 2.8)),
61  C2_(dimensionedScalar::lookupOrDefault("C2", dict, inv(dimArea), 1.83e9))
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
67 void
70 (
71  volScalarField& coalescenceRate,
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 
80  coalescenceRate +=
81  C1_*(pow(fi.x(), 2.0/3.0) + pow(fj.x(), 2.0/3.0))
82  *sqrt(pow(fi.x(), 2.0/9.0) + pow(fj.x(), 2.0/9.0))
84  *exp
85  (
86  - C2_*continuousPhase.thermo().mu()*continuousPhase.rho()
89  /pow3(1 + popBal_.alphas())
90  *pow4(cbrt(fi.x())*cbrt(fj.x())/(cbrt(fi.x()) + cbrt(fj.x())))
91  );
92 }
93 
94 
95 // ************************************************************************* //
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
const dimensionSet dimArea
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
const phaseCompressible::momentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const dimensionSet dimless
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dimensionedScalar exp(const dimensionedScalar &ds)
const phaseModel & continuousPhase() const
Return continuous phase.
dimensionedScalar cbrt(const dimensionedScalar &ds)
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const volScalarField & alphas() const
Return total void of phases belonging to this populationBalance.
dimensionedScalar pow4(const dimensionedScalar &ds)
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
CoulaloglouTavlaridesCoalescence(const populationBalanceModel &popBal, const dictionary &dict)
Namespace for OpenFOAM.