BrownianCollisions.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2019-2025 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "BrownianCollisions.H"
28 #include "fundamentalConstants.H"
29 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace coalescenceModels
37 {
40  (
44  );
45 }
46 }
47 }
48 
51 
52 
53 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
54 
57 (
58  const populationBalanceModel& popBal,
59  const dictionary& dict
60 )
61 :
62  coalescenceModel(popBal, dict),
63  A1_(dict.lookupOrDefault<scalar>("A1", 2.514)),
64  A2_(dict.lookupOrDefault<scalar>("A2", 0.8)),
65  A3_(dict.lookupOrDefault<scalar>("A3", 0.55)),
66  sigma_("sigma", dimLength, dict),
67  lambda_
68  (
69  IOobject
70  (
71  "lambda",
72  popBal_.time().name(),
73  popBal_.mesh()
74  ),
75  popBal_.mesh(),
77  (
78  "lambda",
79  dimLength,
80  Zero
81  )
82  )
83 {}
84 
85 
86 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
87 
89 {
91  popBal_.continuousPhase().fluidThermo().p();
92 
93  const volScalarField::Internal& Tc = popBal_.continuousPhase().thermo().T();
94 
95  lambda_ = k*Tc/(sqrt(2.0)*pi*p*sqr(sigma_));
96 }
97 
98 
101 (
102  volScalarField::Internal& coalescenceRate,
103  const label i,
104  const label j
105 )
106 {
107  const sizeGroup& fi = popBal_.sizeGroups()[i];
108  const sizeGroup& fj = popBal_.sizeGroups()[j];
109 
110  tmp<volScalarField> tdi = fi.d();
111  const volScalarField::Internal& di = tdi();
112  tmp<volScalarField> tdj = fj.d();
113  const volScalarField::Internal& dj = tdj();
114 
115  const volScalarField::Internal& Tc = popBal_.continuousPhase().thermo().T();
116 
117  tmp<volScalarField> tmuc(popBal_.continuousPhase().fluidThermo().mu());
118  const volScalarField::Internal& muc = tmuc();
119 
120  const volScalarField::Internal Cci
121  (
122  1 + lambda_/di*(A1_ + A2_*exp(-A3_*di/lambda_))
123  );
124 
125  const volScalarField::Internal Ccj
126  (
127  1 + lambda_/dj*(A1_ + A2_*exp(-A3_*dj/lambda_))
128  );
129 
130  coalescenceRate += 8*k*Tc/(3*muc)*(di + dj)*(Cci/di + Ccj/dj);
131 }
132 
133 
134 // ************************************************************************* //
label k
Macros for easy insertion into run-time selection tables.
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
tmp< DimensionedField< Type, GeoMesh, Field > > T() const
Return the field transpose (only defined for second rank tensors)
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Base class for coalescence models.
Model describing coagulation due to Brownian motion. Utilises collisional diameters and the Cunningha...
virtual void precompute()
Precompute diameter independent expressions.
virtual void addToCoalescenceRate(volScalarField::Internal &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
BrownianCollisions(const populationBalanceModel &popBal, const dictionary &dict)
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:96
const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroup.C:147
A list of keywords followed by any number of values (e.g. words and numbers) or sub-dictionaries.
Definition: dictionary.H:162
A class for managing temporary objects.
Definition: tmp.H:55
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
Fundamental dimensioned constants.
const dimensionedScalar k
Boltzmann constant.
addToRunTimeSelectionTable(coalescenceModel, AdachiStuartFokkink, dictionary)
defineTypeNameAndDebug(AdachiStuartFokkink, 0)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
dimensionedScalar exp(const dimensionedScalar &ds)
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:59
const dimensionSet dimLength
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
void sqrt(LagrangianPatchField< scalar > &f, const LagrangianPatchField< scalar > &f1)
dictionary dict
volScalarField & p