All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
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 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"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
33 namespace diameterModels
34 {
35 namespace coalescenceModels
36 {
37  defineTypeNameAndDebug(AdachiStuartFokkink, 0);
39  (
40  coalescenceModel,
41  AdachiStuartFokkink,
42  dictionary
43  );
44 }
45 }
46 }
47 
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const populationBalanceModel& popBal,
57  const dictionary& dict
58 )
59 :
60  coalescenceModel(popBal, dict)
61 {}
62 
63 
64 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
65 
68 (
69  volScalarField& coalescenceRate,
70  const label i,
71  const label j
72 )
73 {
74  const sizeGroup& fi = popBal_.sizeGroups()[i];
75  const sizeGroup& fj = popBal_.sizeGroups()[j];
76 
77  coalescenceRate +=
78  (4.0/3.0)
79  *sqrt
80  (
83  )
84  *pow3(fi.d() + fj.d());
85 }
86 
87 
88 // ************************************************************************* //
dictionary dict
FvWallInfoData< WallInfo, label > label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
const phaseCompressible::momentumTransportModel & continuousTurbulence() const
Return reference to turbulence model of the continuous phase.
dimensionedScalar sqrt(const dimensionedScalar &ds)
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
AdachiStuartFokkink(const populationBalanceModel &popBal, const dictionary &dict)
const phaseModel & continuousPhase() const
Return continuous phase.
const rhoThermo & thermo() const
Return const-access to phase rhoThermo.
Definition: phaseModel.H:121
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: breakupModel.H:62
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the size groups belonging to this populationBalance.
dimensionedScalar pow3(const dimensionedScalar &ds)
virtual tmp< volScalarField > nu() const
Kinematic viscosity of mixture [m^2/s].
Definition: fluidThermo.C:73
virtual tmp< volScalarField > epsilon() const =0
Return the turbulence kinetic energy dissipation rate.
Namespace for OpenFOAM.