Loading...
Searching...
No Matches
PrinceBlanch.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
29#include "PrinceBlanch.H"
32#include "phaseCompressibleTurbulenceModel.H"
33
34// * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
36namespace Foam
37{
38namespace diameterModels
39{
40namespace coalescenceModels
41{
44 (
48 );
49}
51}
52
54
55// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
56
59(
60 const populationBalanceModel& popBal,
61 const dictionary& dict
62)
63:
64 coalescenceModel(popBal, dict),
65 C1_(dimensionedScalar::getOrDefault("C1", dict, dimless, 0.356)),
66 h0_
67 (
68 dimensionedScalar::getOrDefault
69 (
70 "h0",
71 dict,
73 1e-4
74 )
75 ),
76 hf_
77 (
78 dimensionedScalar::getOrDefault
79 (
80 "hf",
81 dict,
83 1e-8
84 )
85 ),
86 turbulence_(dict.lookup("turbulence")),
87 buoyancy_(dict.lookup("buoyancy")),
88 laminarShear_(dict.lookup("laminarShear"))
89{}
90
91
92// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93
96(
97 volScalarField& coalescenceRate,
98 const label i,
99 const label j
100)
101{
102 const phaseModel& continuousPhase = popBal_.continuousPhase();
103 const sizeGroup& fi = popBal_.sizeGroups()[i];
104 const sizeGroup& fj = popBal_.sizeGroups()[j];
107
108 const dimensionedScalar rij(1.0/(1.0/fi.d() + 1.0/fj.d()));
109
110 const volScalarField collisionEfficiency
111 (
112 exp
113 (
114 - sqrt
115 (
116 pow3(rij)*continuousPhase.rho()
118 )
119 *log(h0_/hf_)
120 *cbrt(popBal_.continuousTurbulence().epsilon())/pow(rij, 2.0/3.0)
121 )
122 );
123
124 if (turbulence_)
125 {
126 coalescenceRate +=
127 (
128 C1_*pi*sqr(fi.d() + fj.d())
130 *sqrt(pow(fi.d(), 2.0/3.0) + pow(fj.d(), 2.0/3.0))
131 )
132 *collisionEfficiency;
133 }
134
135 if (buoyancy_)
136 {
137 const dimensionedScalar Sij(pi/4.0*sqr(fi.d() + fj.d()));
138
139 coalescenceRate +=
140 (
141 Sij
142 *mag
143 (
144 sqrt
145 (
147 /(continuousPhase.rho()*fi.d()) + 0.505*mag(g)*fi.d()
148 )
149 - sqrt
150 (
152 /(continuousPhase.rho()*fj.d()) + 0.505*mag(g)*fj.d()
153 )
154 )
155 )
156 *collisionEfficiency;
157 }
158
159 if (laminarShear_)
160 {
162 << "Laminar shear collision contribution not implemented for "
163 << this->type() << " coalescence model."
164 << exit(FatalError);
165 }
166}
167
168
169// ************************************************************************* //
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.
const uniformDimensionedVectorField & g
Base class for coalescence models.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
coalescenceModel(const populationBalanceModel &popBal, const dictionary &dict)
Model of Prince and Blanch (1990). The coalescence rate is calculated by.
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
PrinceBlanch(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 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.
const fvMesh & mesh() const
Return reference to the mesh.
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
const Time & time() const
Return the top-level database.
Definition fvMesh.H:360
const Type & lookupObject(const word &name, const bool recursive=false) const
Lookup and return const reference to the object of the given Type. Fatal if not found or the wrong ty...
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
Lookup type of boundary radiation properties.
Definition lookup.H:60
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
#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)
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
dimensionedScalar pow3(const dimensionedScalar &ds)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
fileName::Type type(const fileName &name, const bool followLink=true)
Return the file type: DIRECTORY or FILE, normally following symbolic links.
Definition POSIX.C:801
dimensionedScalar log(const dimensionedScalar &ds)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
UniformDimensionedField< vector > uniformDimensionedVectorField
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
volScalarField & e