populationBalanceSystem.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) 2017-2025 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 
27 #include "fvmSup.H"
28 
29 // * * * * * * * * * * * * * * Static Data Members * * * * * * * * * * * * * //
30 
31 namespace Foam
32 {
34 }
35 
36 
37 // * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * //
38 
39 void Foam::populationBalanceSystem::addDmdts
40 (
41  const diameterModels::populationBalanceModel::dmdtfTable& dmdtfs,
42  PtrList<volScalarField::Internal>& dmdts
43 ) const
44 {
45  forAll(populationBalances_, popBali)
46  {
48  (
50  dmdtfs,
51  dmdtfIter
52  )
53  {
54  const phaseInterface interface(fluid_, dmdtfIter.key());
55 
56  addField(interface.phase1(), "dmdt", *dmdtfIter(), dmdts);
57  addField(interface.phase2(), "dmdt", - *dmdtfIter(), dmdts);
58  }
59  }
60 }
61 
62 
63 void Foam::populationBalanceSystem::addDmdtUfs
64 (
66  HashPtrTable<fvVectorMatrix>& eqns
67 ) const
68 {
70  (
72  dmdtfs,
73  dmdtfIter
74  )
75  {
76  const phaseInterface interface(fluid_, dmdtfIter.key());
77 
78  const volScalarField::Internal& dmdtf = *dmdtfIter();
79  const volScalarField::Internal dmdtf21(posPart(dmdtf));
80  const volScalarField::Internal dmdtf12(negPart(dmdtf));
81 
82  const phaseModel& phase1 = fluid_.phases()[interface.phase1().name()];
83  const phaseModel& phase2 = fluid_.phases()[interface.phase2().name()];
84 
85  if (!phase1.stationary())
86  {
87  *eqns[phase1.name()] +=
88  dmdtf21*phase2.U()()() + fvm::Sp(dmdtf12, phase1.URef());
89  }
90 
91  if (!phase2.stationary())
92  {
93  *eqns[phase2.name()] -=
94  dmdtf12*phase1.U()()() + fvm::Sp(dmdtf21, phase2.URef());
95  }
96  }
97 }
98 
99 
100 void Foam::populationBalanceSystem::addDmdtHefs
101 (
103  HashPtrTable<fvScalarMatrix>& eqns
104 ) const
105 {
106  // Loop the pairs
108  (
110  dmdtfs,
111  dmdtfIter
112  )
113  {
114  const phaseInterface interface(fluid_, dmdtfIter.key());
115 
116  const volScalarField::Internal& dmdtf = *dmdtfIter();
117  const volScalarField::Internal dmdtf21(posPart(dmdtf));
118  const volScalarField::Internal dmdtf12(negPart(dmdtf));
119 
120  const phaseModel& phase1 = interface.phase1();
121  const phaseModel& phase2 = interface.phase2();
122  const rhoFluidThermo& thermo1 = phase1.fluidThermo();
123  const rhoFluidThermo& thermo2 = phase2.fluidThermo();
124  const volScalarField& he1 = thermo1.he();
125  const volScalarField& he2 = thermo2.he();
126  const volScalarField::Internal hs1(thermo1.hs());
127  const volScalarField::Internal hs2(thermo2.hs());
128  const volScalarField::Internal K1(phase1.K());
129  const volScalarField::Internal K2(phase2.K());
130 
131  // Transfer of sensible enthalpy within the phases
132  *eqns[phase1.name()] +=
133  dmdtf*hs1 + fvm::Sp(dmdtf12, he1) - dmdtf12*he1;
134  *eqns[phase2.name()] -=
135  dmdtf*hs2 + fvm::Sp(dmdtf21, he2) - dmdtf21*he2;
136 
137  // Transfer of sensible enthalpy between the phases
138  *eqns[phase1.name()] += dmdtf21*(hs2 - hs1);
139  *eqns[phase2.name()] -= dmdtf12*(hs1 - hs2);
140 
141  // Transfer of kinetic energy
142  *eqns[phase1.name()] += dmdtf21*K2 + dmdtf12*K1;
143  *eqns[phase2.name()] -= dmdtf12*K1 + dmdtf21*K2;
144  }
145 }
146 
147 
148 void Foam::populationBalanceSystem::addDmdtYfs
149 (
151  HashPtrTable<fvScalarMatrix>& eqns
152 ) const
153 {
155  (
157  dmdtfs,
158  dmdtfIter
159  )
160  {
161  const phaseInterface interface(fluid_, dmdtfIter.key());
162 
163  const volScalarField::Internal& dmdtf = *dmdtfIter();
164  const volScalarField::Internal dmdtf12(negPart(dmdtf));
165  const volScalarField::Internal dmdtf21(posPart(dmdtf));
166 
167  const phaseModel& phase1 = interface.phase1();
168  const phaseModel& phase2 = interface.phase2();
169 
170  forAll(phase1.Y(), Yi1)
171  {
172  const volScalarField& Y1 = phase1.Y()[Yi1];
173  const volScalarField& Y2 = phase2.Y(Y1.member());
174 
175  *eqns[Y1.name()] += dmdtf21*Y2 + fvm::Sp(dmdtf12, Y1);
176  *eqns[Y2.name()] -= dmdtf12*Y1 + fvm::Sp(dmdtf21, Y2);
177  }
178  }
179 }
180 
181 
182 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
183 
185 (
186  const phaseSystem& fluid
187 )
188 :
189  fluid_(fluid),
190  populationBalances_()
191 {
192  if (fluid_.found("populationBalances"))
193  {
195  (
196  fluid_.lookup("populationBalances"),
198  );
199 
200  populationBalances_.transfer(popBals);
201  }
202 }
203 
204 
205 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
206 
208 {}
209 
210 
211 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
212 
215 {
216  PtrList<volScalarField::Internal> dmdts(fluid_.phases().size());
217 
218  forAll(populationBalances_, popBali)
219  {
220  addDmdts(populationBalances_[popBali].dmdtfs(), dmdts);
221  addDmdts(populationBalances_[popBali].expansionDmdtfs(), dmdts);
222  addDmdts(populationBalances_[popBali].modelSourceDmdtfs(), dmdts);
223  }
224 
225  return dmdts;
226 }
227 
228 
231 {
233  (
235  );
236  HashPtrTable<fvVectorMatrix>& eqns = eqnsPtr();
237 
238  forAll(fluid_.movingPhases(), movingPhasei)
239  {
240  const phaseModel& phase = fluid_.movingPhases()[movingPhasei];
241 
242  eqns.insert
243  (
244  phase.name(),
246  );
247  }
248 
249  forAll(populationBalances_, popBali)
250  {
251  addDmdtUfs(populationBalances_[popBali].dmdtfs(), eqns);
252  addDmdtUfs(populationBalances_[popBali].expansionDmdtfs(), eqns);
253  addDmdtUfs(populationBalances_[popBali].modelSourceDmdtfs(), eqns);
254  }
255 
256  return eqnsPtr;
257 }
258 
259 
262 {
263  return momentumTransfer();
264 }
265 
266 
269 {
271  (
273  );
274  HashPtrTable<fvScalarMatrix>& eqns = eqnsPtr();
275 
276  forAll(fluid_.phases(), phasei)
277  {
278  const phaseModel& phase = fluid_.phases()[phasei];
279 
280  eqns.insert
281  (
282  phase.name(),
283  new fvScalarMatrix(phase.thermo().he(), dimEnergy/dimTime)
284  );
285  }
286 
287  forAll(populationBalances_, popBali)
288  {
289  addDmdtHefs(populationBalances_[popBali].dmdtfs(), eqns);
290  addDmdtHefs(populationBalances_[popBali].expansionDmdtfs(), eqns);
291  addDmdtHefs(populationBalances_[popBali].modelSourceDmdtfs(), eqns);
292  }
293 
294  return eqnsPtr;
295 }
296 
297 
300 {
302  (
304  );
305  HashPtrTable<fvScalarMatrix>& eqns = eqnsPtr();
306 
307  forAll(fluid_.multicomponentPhases(), multicomponentPhasei)
308  {
309  const phaseModel& phase =
310  fluid_.multicomponentPhases()[multicomponentPhasei];
311 
312  const UPtrList<volScalarField>& Y = phase.Y();
313 
314  forAll(Y, i)
315  {
316  eqns.insert
317  (
318  Y[i].name(),
320  );
321  }
322  }
323 
324  forAll(populationBalances_, popBali)
325  {
326  addDmdtYfs(populationBalances_[popBali].dmdtfs(), eqns);
327  addDmdtYfs(populationBalances_[popBali].expansionDmdtfs(), eqns);
328  addDmdtYfs(populationBalances_[popBali].modelSourceDmdtfs(), eqns);
329  }
330 
331  return eqnsPtr;
332 }
333 
334 
336 {
337  forAll(populationBalances_, i)
338  {
339  populationBalances_[i].solve();
340  }
341 }
342 
343 
345 {
346  forAll(populationBalances_, i)
347  {
348  populationBalances_[i].correct();
349  }
350 }
351 
352 
353 // ************************************************************************* //
#define K1
Definition: SHA1.C:167
#define K2
Definition: SHA1.C:168
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:476
DimensionedField< Type, GeoMesh, PrimitiveField > Internal
Type of the internal field from which this GeometricField is derived.
A HashTable specialisation for hashing pointers.
Definition: HashPtrTable.H:67
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: UPtrList.H:66
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
virtual const volScalarField & he() const =0
Enthalpy/Internal energy [J/kg].
Return a pointer to a new populationBalanceModel object created on.
HashPtrTable< volScalarField::Internal, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
Table of interfacial mass transfer rates.
ITstream & lookup(const word &, bool recursive=false, bool patternMatch=true) const
Find and return an entry data stream.
Definition: dictionary.C:740
bool found(const word &, bool recursive=false, bool patternMatch=true) const
Search dictionary for given keyword.
Definition: dictionary.C:539
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvMatrix.H:118
virtual const rhoThermo & thermo() const =0
Return the thermophysical model.
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
virtual tmp< volVectorField > U() const =0
Return the velocity.
const word & name() const
Return the name of this phase.
Definition: phaseModel.C:145
Class to represent a system of phases.
Definition: phaseSystem.H:74
Class which provides population balance functionality. Holds a number of population balances and prov...
void correct()
Correct the population balances.
autoPtr< HashPtrTable< fvVectorMatrix > > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
autoPtr< HashPtrTable< fvScalarMatrix > > specieTransfer() const
Return the specie transfer matrices.
autoPtr< HashPtrTable< fvScalarMatrix > > heatTransfer() const
Return the heat transfer matrices.
populationBalanceSystem(const phaseSystem &)
Construct from a phase system.
PtrList< volScalarField::Internal > dmdts() const
Return the mass transfer rates for each phase.
autoPtr< HashPtrTable< fvVectorMatrix > > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
void solve()
Solve all population balance equations.
virtual ~populationBalanceSystem()
Destructor.
Calculate the matrix for implicit and explicit sources.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
Namespace for OpenFOAM.
const dimensionSet dimEnergy
void addField(const label phasei, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
Definition: phaseSystemI.H:265
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
const dimensionSet dimTime
dimensionedScalar negPart(const dimensionedScalar &ds)
VolField< scalar > volScalarField
Definition: volFieldsFwd.H:62
defineTypeNameAndDebug(combustionModel, 0)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
const dimensionSet dimVelocity
dimensionedScalar posPart(const dimensionedScalar &ds)
PtrList< volScalarField > & Y