AdachiStuartFokkink.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) 2022-2023 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 "AdachiStuartFokkink.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 
50 
51 
52 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
53 
56 (
57  const populationBalanceModel& popBal,
58  const dictionary& dict
59 )
60 :
61  coalescenceModel(popBal, dict)
62 {}
63 
64 
65 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
66 
69 (
70  volScalarField& coalescenceRate,
71  const label i,
72  const label j
73 )
74 {
75  const sizeGroup& fi = popBal_.sizeGroups()[i];
76  const sizeGroup& fj = popBal_.sizeGroups()[j];
77 
78  coalescenceRate +=
79  (4.0/3.0)
80  *sqrt
81  (
82  0.3*pi*popBal_.continuousTurbulence().epsilon()
83  /popBal_.continuousPhase().fluidThermo().nu()
84  )
85  *pow3(fi.d() + fj.d());
86 }
87 
88 
89 // ************************************************************************* //
Macros for easy insertion into run-time selection tables.
Generic GeometricField class.
Base class for coalescence models.
Model describing aggregation of solid particles in turbulent flows. Applicable when particles are sma...
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
AdachiStuartFokkink(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:102
const tmp< volScalarField > d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroup.C:142
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
addToRunTimeSelectionTable(coalescenceModel, AdachiStuartFokkink, dictionary)
defineTypeNameAndDebug(AdachiStuartFokkink, 0)
Namespace for OpenFOAM.
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
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar sqrt(const dimensionedScalar &ds)
dictionary dict