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 -------------------------------------------------------------------------------
10 License
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 
30 #include "phaseCompressibleTurbulenceModel.H"
31 
32 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
33 
34 namespace Foam
35 {
36 namespace diameterModels
37 {
38 namespace coalescenceModels
39 {
40  defineTypeNameAndDebug(CoulaloglouTavlaridesCoalescence, 0);
42  (
43  coalescenceModel,
44  CoulaloglouTavlaridesCoalescence,
45  dictionary
46  );
47 }
48 }
49 }
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  coalescenceModel(popBal, dict),
62  C1_(dimensionedScalar::lookupOrDefault("C1", dict, dimless, 2.8)),
63  C2_(dimensionedScalar::lookupOrDefault("C2", dict, inv(dimArea), 1.83e9))
64 {}
65 
66 
67 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
68 
69 void
72 (
73  volScalarField& coalescenceRate,
74  const label i,
75  const label j
76 )
77 {
78  const phaseModel& continuousPhase = popBal_.continuousPhase();
79  const sizeGroup& fi = popBal_.sizeGroups()[i];
80  const sizeGroup& fj = popBal_.sizeGroups()[j];
81 
82  coalescenceRate +=
83  C1_*(pow(fi.x(), 2.0/3.0) + pow(fj.x(), 2.0/3.0))
84  *sqrt(pow(fi.x(), 2.0/9.0) + pow(fj.x(), 2.0/9.0))
85  *cbrt(popBal_.continuousTurbulence().epsilon())/(1 + popBal_.alphas())
86  *exp
87  (
88  - C2_*continuousPhase.mu()*continuousPhase.rho()
89  *popBal_.continuousTurbulence().epsilon()
90  /sqr(popBal_.sigmaWithContinuousPhase(fi.phase()))
91  /pow3(1 + popBal_.alphas())
92  *pow4(cbrt(fi.x())*cbrt(fj.x())/(cbrt(fi.x()) + cbrt(fj.x())))
93  );
94 }
95 
96 
97 // ************************************************************************* //
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:57
Foam::dimless
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Dimensionless.
Definition: dimensionSets.H:50
Foam::diameterModels::sizeGroup::phase
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:38
Foam::phaseModel::rho
virtual tmp< volScalarField > rho() const =0
Return the density field.
Foam::diameterModels::coalescenceModels::CoulaloglouTavlaridesCoalescence::addToCoalescenceRate
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition: CoulaloglouTavlaridesCoalescence.C:72
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:56
CoulaloglouTavlaridesCoalescence.H
Foam::pow4
dimensionedScalar pow4(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:100
Foam::dimArea
const dimensionSet dimArea(sqr(dimLength))
Definition: dimensionSets.H:60
Foam::label
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:62
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:121
Foam::diameterModels::populationBalanceModel
Class that solves the univariate population balance equation by means of a class method (also called ...
Definition: populationBalanceModel.H:179
Foam::phaseModel::mu
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
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::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