All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
DahnekeInterpolation.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-2020 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 "DahnekeInterpolation.H"
28 #include "BrownianCollisions.H"
29 #include "ballisticCollisions.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace coalescenceModels
38 {
39  defineTypeNameAndDebug(DahnekeInterpolation, 0);
41  (
42  coalescenceModel,
43  DahnekeInterpolation,
44  dictionary
45  );
46 }
47 }
48 }
49 
50 
51 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
52 
55 (
56  const populationBalanceModel& popBal,
57  const dictionary& dict
58 )
59 :
60  coalescenceModel(popBal, dict),
61  Brownian_(new BrownianCollisions(popBal, dict)),
62  BrownianRate_
63  (
64  IOobject
65  (
66  "BrownianCollisionRate",
67  popBal_.mesh().time().timeName(),
68  popBal_.mesh()
69  ),
70  popBal_.mesh(),
71  dimensionedScalar("BrownianCollisionRate", dimVolume/dimTime, Zero)
72  ),
73  ballistic_(new ballisticCollisions(popBal, dict)),
74  ballisticRate_
75  (
76  IOobject
77  (
78  "ballisticCollisionRate",
79  popBal_.mesh().time().timeName(),
80  popBal_.mesh()
81  ),
82  popBal_.mesh(),
83  dimensionedScalar("ballisticCollisionRate", dimVolume/dimTime, Zero)
84  )
85 {}
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  Brownian_().precompute();
93 }
94 
95 
96 void
99 (
100  volScalarField& coalescenceRate,
101  const label i,
102  const label j
103 )
104 {
105  BrownianRate_ = Zero;
106  ballisticRate_ = Zero;
107 
108  Brownian_().addToCoalescenceRate(BrownianRate_, i, j);
109  ballistic_().addToCoalescenceRate(ballisticRate_, i, j);
110 
111  const volScalarField KnD(BrownianRate_/(2*ballisticRate_));
112 
113  coalescenceRate += BrownianRate_*(1 + KnD)/(1 + 2*KnD + 2*sqr(KnD));
114 }
115 
116 
117 // ************************************************************************* //
virtual void precompute()
Precompute diameter independent expressions.
dictionary dict
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Macros for easy insertion into run-time selection tables.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const dimensionSet dimTime
dynamicFvMesh & mesh
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimVolume
Namespace for OpenFOAM.
virtual void addToCoalescenceRate(volScalarField &coalescenceRate, const label i, const label j)
Add to coalescenceRate.
DahnekeInterpolation(const populationBalanceModel &popBal, const dictionary &dict)