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-2019 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 {
39  defineTypeNameAndDebug(phaseChange, 0);
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 {
59  const phaseSystem& fluid = popBal_.fluid();
60 
61  forAll(pairKeys_, i)
62  {
63  const phasePair& pair = fluid.phasePairs()[pairKeys_[i]];
64 
65  W_.set
66  (
67  i,
68  new volScalarField
69  (
70  IOobject
71  (
72  IOobject::groupName(type() + ":W", pair.name()),
73  popBal_.mesh().time().timeName(),
74  popBal_.mesh()
75  ),
76  popBal_.mesh(),
78  (
79  inv(numberWeighted_ ? dimVolume : dimLength),
80  Zero
81  )
82  )
83  );
84  }
85 }
86 
87 
88 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
89 
91 {
92  const phaseSystem& fluid = popBal_.fluid();
93 
94  forAll(pairKeys_, i)
95  {
96  W_[i] = Zero;
97  }
98 
99  forAll(pairKeys_, k)
100  {
101  if (fluid.phasePairs().found(pairKeys_[k]))
102  {
103  const phasePair& pair = fluid.phasePairs()[pairKeys_[k]];
104 
106  {
107  const velocityGroup& vgj = popBal_.velocityGroups()[j];
108  if (pair.contains(vgj.phase()))
109  {
110  forAll(vgj.sizeGroups(), i)
111  {
112  const sizeGroup& fi = vgj.sizeGroups()[i];
113 
114  W_[k] +=
115  fi*max(fi.phase(), small)
116  /(numberWeighted_ ? fi.x() : fi.d());
117  }
118  }
119  }
120  }
121  }
122 }
123 
124 
126 (
127  volScalarField& driftRate,
128  const label i
129 )
130 {
131  const velocityGroup& vg = popBal_.sizeGroups()[i].VelocityGroup();
132 
133  forAll(pairKeys_, k)
134  {
135  const phasePair& pair =
136  popBal_.fluid().phasePairs()[pairKeys_[k]];
137 
138  if (pair.contains(vg.phase()))
139  {
140  const volScalarField& iDmdt =
142  (
143  IOobject::groupName("iDmdt", pair.name())
144  );
145 
146  const scalar iDmdtSign =
147  vg.phase().name() == pair.first() ? +1 : -1;
148 
149  const sizeGroup& fi = popBal_.sizeGroups()[i];
150 
151  tmp<volScalarField> dDriftRate
152  (
153  iDmdtSign*iDmdt/(fi.phase().rho()*W_[k])
154  );
155 
156  if (!numberWeighted_)
157  {
158  dDriftRate.ref() *= fi.x()/fi.d();
159  }
160 
161  driftRate += dDriftRate;
162  }
163  }
164 }
165 
166 
167 // ************************************************************************* //
#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 correct()
Correct diameter independent expressions.
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: driftModel.H:57
static word timeName(const scalar, const int precision=precision_)
Return time name of given scalar time.
Definition: Time.C:626
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.
phaseChange(const populationBalanceModel &popBal, const dictionary &dict)
Construct from a population balance model and a dictionary.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:52
const phasePairTable & phasePairs() const
Return the phase pairs.
Definition: phaseSystemI.H:105
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
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.
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:50
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 fvMesh & mesh() const
Return reference to the mesh.
Namespace for OpenFOAM.