Loading...
Searching...
No Matches
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-------------------------------------------------------------------------------
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
31#include "phaseCompressibleTurbulenceModel.H"
32
33// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
35namespace Foam
36{
37namespace diameterModels
38{
39namespace coalescenceModels
40{
43 (
47 );
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
70void
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))
87 *exp
88 (
89 - C2_*continuousPhase.mu()*continuousPhase.rho()
92 /pow3(1 + popBal_.alphas())
93 *pow4(cbrt(fi.x())*cbrt(fj.x())/(cbrt(fi.x()) + cbrt(fj.x())))
94 );
95}
96
97
98// ************************************************************************* //
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.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
coalescenceModel(const populationBalanceModel &popBal, const dictionary &dict)
Model of Coulaloglou and Tavlarides (1977). The coalescence rate is calculated by.
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
CoulaloglouTavlaridesCoalescence(const populationBalanceModel &popBal, const dictionary &dict)
Class that solves the univariate population balance equation by means of a class method (also called ...
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
const phaseModel & continuousPhase() const
Return continuous phase.
const volScalarField & alphas() const
Return total void of phases belonging to this populationBalance.
const phaseCompressibleTurbulenceModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
const tmp< volScalarField > sigmaWithContinuousPhase(const phaseModel &dispersedPhase) const
Return the surface tension coefficient between a given dispersed.
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition sizeGroup.H:95
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition sizeGroupI.H:52
const phaseModel & phase() const
Return const-reference to the phase.
Definition sizeGroupI.H:31
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition dictionary.H:133
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition phaseModel.H:56
const dimensionedScalar & rho() const
Definition phaseModel.H:193
virtual tmp< volScalarField > mu() const =0
Return the laminar dynamic viscosity.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
Namespace for OpenFOAM.
dimensionedScalar exp(const dimensionedScalar &ds)
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
const dimensionSet dimArea(sqr(dimLength))
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar pow4(const dimensionedScalar &ds)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
dictionary dict