Loading...
Searching...
No Matches
Luo.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) 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
29#include "Luo.H"
32#include "phaseCompressibleTurbulenceModel.H"
33#include "virtualMassModel.H"
34#include "phaseSystem.H"
35
36// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
38namespace Foam
39{
40namespace diameterModels
41{
42namespace coalescenceModels
43{
46 (
48 Luo,
50 );
51}
53}
54
56
57// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
58
60Luo
61(
62 const populationBalanceModel& popBal,
63 const dictionary& dict
64)
65:
66 coalescenceModel(popBal, dict),
67 beta_(dimensionedScalar::getOrDefault("beta", dict, dimless, 2.05)),
68 C1_(dimensionedScalar::getOrDefault("C1", dict, dimless, 1.0))
69{}
70
71
72// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
73
76(
77 volScalarField& coalescenceRate,
78 const label i,
79 const label j
80)
81{
82 const sizeGroup& fi = popBal_.sizeGroups()[i];
83 const sizeGroup& fj = popBal_.sizeGroups()[j];
84 const phaseModel& continuousPhase = popBal_.continuousPhase();
85
86 if
87 (
89 (
90 fi.phase(),
92 )
93 )
94 {
95 const virtualMassModel& vm =
97 (
98 fi.phase(),
100 );
101
102 const dimensionedScalar xi = fi.d()/fj.d();
103
104 const volScalarField uij
105 (
106 sqrt(beta_)
108 *sqrt(1.0 + pow(xi, -2.0/3.0))
109 );
110
111 coalescenceRate +=
112 pi/4.0*sqr(fi.d() + fj.d())*uij
113 *exp
114 (
115 - C1_
116 *sqrt(0.75*(1.0 + sqr(xi))*(1.0 + pow3(xi)))
117 /(
118 sqrt(fi.phase().rho()/continuousPhase.rho()
119 + vm.Cvm())*pow3(1.0 + xi)
120 )
121 *sqrt
122 (
123 continuousPhase.rho()*fi.d()*sqr(uij)
125 )
126 );
127 }
128 else
129 {
131 << "A virtual mass model for " << fi.phase().name() << " in "
132 << popBal_.continuousPhase().name() << " is not specified. This is "
133 << "required by the Luo coalescence model." << exit(FatalError);
134 }
135}
136
137
138// ************************************************************************* //
constexpr scalar pi(M_PI)
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 Luo (1993). The coalescence rate is calculated by.
Definition Luo.H:167
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Definition Luo.C:69
Luo(const populationBalanceModel &popBal, const dictionary &dict)
Definition Luo.C:54
Class that solves the univariate population balance equation by means of a class method (also called ...
const phaseSystem & fluid() const
Return reference to the phaseSystem.
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
const phaseModel & continuousPhase() const
Return continuous phase.
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 & d() const
Return representative diameter of the sizeGroup.
Definition sizeGroupI.H:45
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
const word & name() const
Definition phaseModel.H:166
const modelType & lookupSubModel(const phasePair &key) const
Return a sub model between a phase pair.
bool foundSubModel(const phasePair &key) const
Check availability of a sub model for a given phase pair.
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
virtual tmp< volScalarField > Cvm() const =0
Return the virtual mass coefficient.
#define defineTypeNameAndDebug(Type, DebugSwitch)
Define the typeName and debug information.
Definition className.H:142
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition error.H:600
constexpr scalar pi(M_PI)
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
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
error FatalError
Error stream (stdout output on all processes), with additional 'FOAM FATAL ERROR' header text and sta...
dimensionedScalar cbrt(const dimensionedScalar &ds)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition errorManip.H:125
dictionary dict