CoulaloglouTavlaridesCoalescence.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) 2018-2019 OpenFOAM Foundation
9  Copyright (C) 2020 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
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 
31 #include "phaseCompressibleTurbulenceModel.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace diameterModels
38 {
39 namespace coalescenceModels
40 {
41  defineTypeNameAndDebug(CoulaloglouTavlaridesCoalescence, 0);
43  (
44  coalescenceModel,
45  CoulaloglouTavlaridesCoalescence,
46  dictionary
47  );
48 }
49 }
50 }
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
57 (
58  const populationBalanceModel& popBal,
59  const dictionary& dict
60 )
61 :
62  coalescenceModel(popBal, dict),
63  C1_(dimensionedScalar::getOrDefault("C1", dict, dimless, 2.8)),
64  C2_(dimensionedScalar::getOrDefault("C2", dict, inv(dimArea), 1.83e9))
65 {}
66 
67 
68 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
69 
70 void
73 (
74  volScalarField& coalescenceRate,
75  const label i,
76  const label j
77 )
78 {
79  const phaseModel& continuousPhase = popBal_.continuousPhase();
80  const sizeGroup& fi = popBal_.sizeGroups()[i];
81  const sizeGroup& fj = popBal_.sizeGroups()[j];
82 
83  coalescenceRate +=
84  C1_*(pow(fi.x(), 2.0/3.0) + pow(fj.x(), 2.0/3.0))
85  *sqrt(pow(fi.x(), 2.0/9.0) + pow(fj.x(), 2.0/9.0))
86  *cbrt(popBal_.continuousTurbulence().epsilon())/(1 + popBal_.alphas())
87  *exp
88  (
89  - C2_*continuousPhase.mu()*continuousPhase.rho()
90  *popBal_.continuousTurbulence().epsilon()
91  /sqr(popBal_.sigmaWithContinuousPhase(fi.phase()))
92  /pow3(1 + popBal_.alphas())
93  *pow4(cbrt(fi.x())*cbrt(fj.x())/(cbrt(fi.x()) + cbrt(fj.x())))
94  );
95 }
96 
97 
98 // ************************************************************************* //
Foam::phaseModel::mu
virtual tmp< volScalarField > mu() const
Return the mixture dymanic viscosity.
Definition: phaseModel.C:297
Foam::diameterModels::coalescenceModel
Base class for coalescence models.
Definition: coalescenceModel.H:52
Foam::phaseModel
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:54
Foam::diameterModels::sizeGroup::phase
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:38
Foam::diameterModels::coalescenceModels::CoulaloglouTavlaridesCoalescence::addToCoalescenceRate
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition: CoulaloglouTavlaridesCoalescence.C:73
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::diameterModels::coalescenceModels::CoulaloglouTavlaridesCoalescence::CoulaloglouTavlaridesCoalescence
CoulaloglouTavlaridesCoalescence(const populationBalanceModel &popBal, const dictionary &dict)
Definition: CoulaloglouTavlaridesCoalescence.C:57
CoulaloglouTavlaridesCoalescence.H
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:59
Foam::diameterModels::coalescenceModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(coalescenceModel, constantCoalescence, dictionary)
Foam::pow3
dimensionedScalar pow3(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:89
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::diameterModels::sizeGroup
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:96
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::diameterModels::populationBalanceModel
Class that solves the univariate population balance equation by means of a class method (also called ...
Definition: populationBalanceModel.H:179
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:144
Foam::diameterModels::coalescenceModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constantCoalescence, 0)
Foam::phaseModel::rho
const dimensionedScalar & rho() const
Definition: phaseModel.H:167
Foam::cbrt
dimensionedScalar cbrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:155
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::diameterModels::sizeGroup::x
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:59
Foam::dimless
const dimensionSet dimless
Dimensionless.
Definition: dimensionSets.C:189