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-2023 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 
30 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
31 
32 namespace Foam
33 {
34 namespace diameterModels
35 {
36 namespace driftModels
37 {
40 }
41 }
42 }
43 
44 
45 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
46 
48 (
49  const populationBalanceModel& popBal,
50  const dictionary& dict
51 )
52 :
53  driftModel(popBal, dict),
54  interfaces_
55  (
56  dict.lookup("interfaces"),
57  phaseInterface::iNew(popBal_.fluid())
58  ),
59  numberWeighted_(dict.lookupOrDefault<Switch>("numberWeighted", false)),
60  W_(interfaces_.size()),
61  dmdtfName_(dict.lookup("dmdtf")),
62  specieName_(dict.lookupOrDefault("specie", word()))
63 {
64  forAll(interfaces_, i)
65  {
66  W_.set
67  (
68  i,
69  new volScalarField
70  (
71  IOobject
72  (
73  IOobject::groupName(typedName("W"), interfaces_[i].name()),
74  popBal_.mesh().time().name(),
75  popBal_.mesh()
76  ),
77  popBal_.mesh(),
79  (
80  inv(numberWeighted_ ? dimVolume : dimLength),
81  Zero
82  )
83  )
84  );
85  }
86 }
87 
88 
89 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
90 
92 {
93  forAll(interfaces_, i)
94  {
95  W_[i] = Zero;
96  }
97 
98  forAll(interfaces_, k)
99  {
100  forAllConstIter(phaseInterface, interfaces_[k], iter)
101  {
102  const phaseModel& phase = iter();
103 
104  if (!isA<velocityGroup>(phase.diameter())) continue;
105 
106  const velocityGroup& velGroup =
107  refCast<const velocityGroup>(phase.diameter());
108 
109  forAll(velGroup.sizeGroups(), i)
110  {
111  const sizeGroup& fi = velGroup.sizeGroups()[i];
112 
113  if (numberWeighted_)
114  {
115  W_[k] += fi*max(fi.phase(), small)/fi.x();
116  }
117  else
118  {
119  W_[k] += fi*max(fi.phase(), small)/fi.x()*fi.a();
120  }
121  }
122  }
123  }
124 }
125 
126 
128 (
129  volScalarField& driftRate,
130  const label i
131 )
132 {
133  const velocityGroup& velGrp = popBal_.sizeGroups()[i].group();
134 
135  forAll(interfaces_, k)
136  {
137  if (interfaces_[k].contains(velGrp.phase()))
138  {
139  const volScalarField& dmidtf =
140  popBal_.mesh().lookupObject<volScalarField>
141  (
143  (
145  (
146  dmdtfName_,
147  specieName_
148  ),
149  interfaces_[k].name()
150  )
151  );
152 
153  const scalar dmidtfSign =
154  interfaces_[k].index(velGrp.phase()) == 0 ? +1 : -1;
155 
156  const sizeGroup& fi = popBal_.sizeGroups()[i];
157 
158  tmp<volScalarField> dDriftRate
159  (
160  dmidtfSign*dmidtf/(fi.phase().rho()*W_[k])
161  );
162 
163  if (!numberWeighted_)
164  {
165  dDriftRate.ref() *= fi.a();
166  }
167 
168  driftRate += velGrp.phase()*dDriftRate;
169  }
170  }
171 }
172 
173 
174 // ************************************************************************* //
label k
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Macros for easy insertion into run-time selection tables.
const Mesh & mesh() const
Return mesh.
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word groupName(Name name, const word &group)
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
const phaseModel & phase() const
Return the phase.
Class used for the read-construction of.
Definition: driftModel.H:90
Base class for drift models.
Definition: driftModel.H:54
const populationBalanceModel & popBal_
Reference to the populationBalanceModel.
Definition: driftModel.H:60
Drift induced by phase change. By default phase change mass flux is distributed between sizeGroups of...
Definition: phaseChange.H:57
virtual void precompute()
Precompute diameter independent expressions.
Definition: phaseChange.C:91
phaseChange(const populationBalanceModel &popBal, const dictionary &dict)
Construct from a population balance model and a dictionary.
Definition: phaseChange.C:48
virtual void addToDriftRate(volScalarField &driftRate, const label i)
Add to driftRate.
Definition: phaseChange.C:128
Model for tracking the evolution of a dispersed phase size distribution due to coalescence (synonymou...
const fvMesh & mesh() const
Return reference to the mesh.
Single size class fraction field representing a fixed particle volume as defined by the user through ...
Definition: sizeGroup.H:102
const dimensionedScalar & x() const
Return representative volume of the sizeGroup.
Definition: sizeGroupI.H:58
const tmp< volScalarField > a() const
Return representative surface area of the sizeGroup.
Definition: sizeGroup.C:135
const phaseModel & phase() const
Return const-reference to the phase.
Definition: sizeGroupI.H:37
Computes the Sauter mean diameter based on a user specified size distribution, defined in terms of si...
Definition: velocityGroup.H:87
const PtrList< sizeGroup > & sizeGroups() const
Return sizeGroups belonging to this velocityGroup.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Definition: dictionary.H:162
const word & name() const
Return const reference to name.
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:418
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const diameterModel & diameter() const
Return a reference to the diameterModel of the phase.
Definition: phaseModel.C:187
virtual const volScalarField & rho() const =0
Return the density field.
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
A class for handling words, derived from string.
Definition: word.H:62
addToRunTimeSelectionTable(driftModel, constantDrift, dictionary)
defineTypeNameAndDebug(constantDrift, 0)
Namespace for OpenFOAM.
static const zero Zero
Definition: zero.H:97
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
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
const dimensionSet dimLength
const dimensionSet dimVolume
word typedName(Name name)
Return the name of the object within the given type.
Definition: typeInfo.H:176
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
dimensionedSphericalTensor inv(const dimensionedSphericalTensor &dt)
dictionary dict