phaseChange.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2018-2019 OpenFOAM Foundation
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
28 #include "phaseChange.H"
30 #include "phaseSystem.H"
31 #include "phasePairKey.H"
32 
33 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
34 
35 namespace Foam
36 {
37 namespace diameterModels
38 {
39 namespace driftModels
40 {
41  defineTypeNameAndDebug(phaseChange, 0);
42  addToRunTimeSelectionTable(driftModel, phaseChange, dictionary);
43 }
44 }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
49 
51 (
52  const populationBalanceModel& popBal,
53  const dictionary& dict
54 )
55 :
56  driftModel(popBal, dict),
57  pairKeys_(dict.lookup("pairs")),
58  numberWeighted_(dict.lookupOrDefault<Switch>("numberWeighted", false)),
59  W_(pairKeys_.size())
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] *= 0.0;
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  W_[k] +=
117  fi*max(fi.phase(), SMALL)
118  /(numberWeighted_ ? fi.x() : fi.d());
119  }
120  }
121  }
122  }
123  }
124 }
125 
126 
128 (
129  volScalarField& driftRate,
130  const label i
131 )
132 {
133  const velocityGroup& vg = popBal_.sizeGroups()[i].VelocityGroup();
134 
135  forAll(pairKeys_, k)
136  {
137  const phasePair& pair =
138  popBal_.fluid().phasePairs()[pairKeys_[k]];
139 
140  if (pair.contains(vg.phase()))
141  {
142  const volScalarField& iDmdt =
143  popBal_.mesh().lookupObject<volScalarField>
144  (
145  IOobject::groupName("iDmdt", pair.name())
146  );
147 
148  const scalar iDmdtSign =
149  vg.phase().name() == pair.first() ? +1 : -1;
150 
151  const sizeGroup& fi = popBal_.sizeGroups()[i];
152 
153  tmp<volScalarField> dDriftRate
154  (
155  iDmdtSign*iDmdt/(fi.phase().rho()*W_[k])
156  );
157 
158  if (!numberWeighted_)
159  {
160  dDriftRate.ref() *= fi.x()/fi.d();
161  }
162 
163  driftRate += dDriftRate;
164  }
165  }
166 }
167 
168 
169 // ************************************************************************* //
Foam::diameterModels::populationBalanceModel::fluid
const phaseSystem & fluid() const
Return reference to the phaseSystem.
Definition: populationBalanceModelI.H:48
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:104
Foam::phasePair
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:51
Foam::Switch
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:70
Foam::diameterModels::sizeGroup::phase
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:38
Foam::dimLength
const dimensionSet dimLength(0, 1, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:53
Foam::tmp< volScalarField >
Foam::Zero
static constexpr const zero Zero
Global zero.
Definition: zero.H:128
Foam::phaseModel::rho
virtual tmp< volScalarField > rho() const =0
Return the density field.
Foam::diameterModels::driftModels::phaseChange::correct
virtual void correct()
Correct diameter independent expressions.
Definition: phaseChange.C:92
fluid
twoPhaseSystem & fluid
Definition: setRegionFluidFields.H:3
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:290
Foam::diameterModels::populationBalanceModel::velocityGroups
const UPtrList< velocityGroup > & velocityGroups() const
Return the velocityGroups belonging to this populationBalance.
Definition: populationBalanceModelI.H:76
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:258
Foam::label
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:62
Foam::inv
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
Definition: dimensionedSphericalTensor.C:73
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:43
Foam::diameterModels::velocityGroup
This diameterModel is intended for use with a populationBalanceModel in order to simulate polydispers...
Definition: velocityGroup.H:107
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
Foam::diameterModels::sizeGroup
This class represents a single sizeGroup belonging to a velocityGroup. The main property of a sizeGro...
Definition: sizeGroup.H:96
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:121
Foam::phasePair::contains
bool contains(const phaseModel &phase) const
Return true if this phasePair contains the given phase.
Definition: phasePairI.H:42
Foam::diameterModels::populationBalanceModel
Class that solves the univariate population balance equation by means of a class method (also called ...
Definition: populationBalanceModel.H:179
Foam::diameterModels::driftModels::defineTypeNameAndDebug
defineTypeNameAndDebug(constantDrift, 0)
addToRunTimeSelectionTable.H
Macros for easy insertion into run-time selection tables.
Foam
Namespace for OpenFOAM.
Definition: atmBoundaryLayer.C:33
Foam::phaseModel::name
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:94
Foam::diameterModels::sizeGroup::d
const dimensionedScalar & d() const
Return representative diameter of the sizeGroup.
Definition: sizeGroupI.H:52
Foam::IOobject::groupName
static word groupName(StringType name, const word &group)
Create dot-delimited name.group.
Foam::phasePair::name
virtual word name() const
Pair name.
Definition: phasePair.C:92
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::roots::type
type
Types of root.
Definition: Roots.H:54
Foam::phaseSystem
Class to represent a system of phases and model interfacial transfers between them.
Definition: phaseSystem.H:69
Foam::diameterModels::driftModel
Base class for drift models.
Definition: driftModel.H:52
Foam::diameterModels::velocityGroup::sizeGroups
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
Definition: velocityGroupI.H:52
Foam::diameterModels::driftModel::popBal_
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: driftModel.H:59
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:61
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::diameterModels::sizeGroup::x
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:59
Foam::diameterModels::driftModels::phaseChange::phaseChange
phaseChange(const populationBalanceModel &popBal, const dictionary &dict)
Construct from a population balance model and a dictionary.
Definition: phaseChange.C:51
Foam::diameterModel::phase
const phaseModel & phase() const
Return the phase.
Definition: diameterModel.H:116
Foam::diameterModels::driftModels::addToRunTimeSelectionTable
addToRunTimeSelectionTable(driftModel, constantDrift, dictionary)
phaseChange.H
Foam::diameterModels::driftModels::phaseChange::addToDriftRate
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
Definition: phaseChange.C:128