phaseChange.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) 2018-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 "phaseChange.H"
28 #include "phaseSystem.H"
29 #include "phasePairKey.H"
30 
31 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
32 
33 namespace Foam
34 {
35 namespace diameterModels
36 {
37 namespace driftModels
38 {
40  addToRunTimeSelectionTable(driftModel, phaseChange, dictionary);
41 }
42 }
43 }
44 
45 
46 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
47 
49 (
50  const populationBalanceModel& popBal,
51  const dictionary& dict
52 )
53 :
54  driftModel(popBal, dict),
55  pairKeys_(dict.lookup("pairs")),
56  numberWeighted_(dict.lookupOrDefault<Switch>("numberWeighted", false)),
57  W_(pairKeys_.size()),
58  dmdtfName_(dict.lookup("dmdtf")),
59  specieName_(dict.lookupOrDefault("specie", word()))
60 {
61  const phaseSystem& fluid = popBal_.fluid();
62 
63  forAll(pairKeys_, i)
64  {
65  const phasePair& pair = fluid.phasePairs()[pairKeys_[i]];
66 
67  W_.set
68  (
69  i,
70  new volScalarField
71  (
72  IOobject
73  (
74  IOobject::groupName(type() + ":W", pair.name()),
75  popBal_.mesh().time().timeName(),
76  popBal_.mesh()
77  ),
78  popBal_.mesh(),
80  (
81  inv(numberWeighted_ ? dimVolume : dimLength),
82  Zero
83  )
84  )
85  );
86  }
87 }
88 
89 
90 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
91 
93 {
94  const phaseSystem& fluid = popBal_.fluid();
95 
96  forAll(pairKeys_, i)
97  {
98  W_[i] = Zero;
99  }
100 
101  forAll(pairKeys_, k)
102  {
103  if (fluid.phasePairs().found(pairKeys_[k]))
104  {
105  const phasePair& pair = fluid.phasePairs()[pairKeys_[k]];
106 
108  {
109  const velocityGroup& vgj = popBal_.velocityGroups()[j];
110  if (pair.contains(vgj.phase()))
111  {
112  forAll(vgj.sizeGroups(), i)
113  {
114  const sizeGroup& fi = vgj.sizeGroups()[i];
115 
116  if (numberWeighted_)
117  {
118  W_[k] += fi*max(fi.phase(), small)/fi.x();
119  }
120  else
121  {
122  W_[k] += fi*max(fi.phase(), small)/fi.x()*fi.a();
123  }
124  }
125  }
126  }
127  }
128  }
129 }
130 
131 
133 (
134  volScalarField& driftRate,
135  const label i
136 )
137 {
138  const velocityGroup& vg = popBal_.sizeGroups()[i].VelocityGroup();
139 
140  forAll(pairKeys_, k)
141  {
142  const phasePair& pair =
143  popBal_.fluid().phasePairs()[pairKeys_[k]];
144 
145  if (pair.contains(vg.phase()))
146  {
147  const volScalarField& dmidtf =
149  (
151  (
153  (
154  dmdtfName_,
155  specieName_
156  ),
157  pair.name()
158  )
159  );
160 
161  const scalar dmidtfSign =
162  vg.phase().name() == pair.first() ? +1 : -1;
163 
164  const sizeGroup& fi = popBal_.sizeGroups()[i];
165 
166  tmp<volScalarField> dDriftRate
167  (
168  dmidtfSign*dmidtf/(fi.phase().rho()*W_[k])
169  );
170 
171  if (!numberWeighted_)
172  {
173  dDriftRate.ref() *= fi.a();
174  }
175 
176  driftRate += dDriftRate;
177  }
178  }
179 }
180 
181 
182 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
virtual void precompute()
Precompute diameter independent expressions.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: driftModel.H:57
twoPhaseChangeModel & phaseChange
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:636
label k
Boltzmann constant.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:239
const Type & lookupObject(const word &name) const
Lookup and return the object of the given Type.
Macros for easy insertion into run-time selection tables.
const dimensionSet dimLength
phaseChange(const populationBalanceModel &popBal, const dictionary &dict)
Construct from a population balance model and a dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
const phasePairTable & phasePairs() const
Return the phase pairs.
Definition: phaseSystemI.H:105
const UPtrList< velocityGroup > & velocityGroups() const
Return the velocityGroups belonging to this populationBalance.
stressControl lookup("compactNormalStress") >> compactNormalStress
const phaseSystem & fluid() const
Return reference to the phaseSystem.
static word groupName(Name name, const word &group)
static const zero Zero
Definition: zero.H:97
addToRunTimeSelectionTable(ensightPart, ensightPartCells, istream)
defineTypeNameAndDebug(combustionModel, 0)
const UPtrList< sizeGroup > & sizeGroups() const
Return the sizeGroups belonging to this populationBalance.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
fileType type(const fileName &, const bool checkVariants=true, const bool followLink=true)
Return the file type: directory or file.
Definition: POSIX.C:488
const dimensionSet dimVolume
const fvMesh & mesh() const
Return reference to the mesh.
Namespace for OpenFOAM.