PhaseTransferPhaseSystem.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 
27 #include "phaseTransferModel.H"
28 #include "fvmSup.H"
29 
30 // * * * * * * * * * * * * Private Member Functions * * * * * * * * * * * * //
31 
32 template<class BasePhaseSystem>
35 {
36  autoPtr<phaseSystem::dmdtfTable> totalDmdtfsPtr
37  (
38  new phaseSystem::dmdtfTable
39  );
40  phaseSystem::dmdtfTable& totalDmdtfs = totalDmdtfsPtr();
41 
43  (
44  phaseTransferModelTable,
45  phaseTransferModels_,
46  phaseTransferModelIter
47  )
48  {
49  const phaseInterface& interface = phaseTransferModelIter()->interface();
50 
51  totalDmdtfs.insert(interface, phaseSystem::dmdtf(interface).ptr());
52 
53  if (phaseTransferModelIter()->mixture())
54  {
55  *totalDmdtfs[interface] += *dmdtfs_[interface];
56  }
57 
59  (
60  HashPtrTable<volScalarField>,
61  *dmidtfs_[interface],
62  dmidtfIter
63  )
64  {
65  *totalDmdtfs[interface] += *dmidtfIter();
66  }
67  }
68 
69  return totalDmdtfsPtr;
70 }
71 
72 
73 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
74 
75 template<class BasePhaseSystem>
77 (
78  const phaseSystem::dmdtfTable& dmdtfs,
80 ) const
81 {
82  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
83  {
84  const phaseInterface interface(*this, dmdtfIter.key());
85 
86  const volScalarField& dmdtf = *dmdtfIter();
87  const volScalarField dmdtf12(negPart(dmdtf));
88  const volScalarField dmdtf21(posPart(dmdtf));
89 
90  const phaseModel& phase1 = interface.phase1();
91  const phaseModel& phase2 = interface.phase2();
92 
93  forAll(phase1.Y(), Yi1)
94  {
95  const volScalarField& Y1 = phase1.Y()[Yi1];
96  const volScalarField& Y2 = phase2.Y(Y1.member());
97 
98  *eqns[Y1.name()] += dmdtf21*Y2 + fvm::Sp(dmdtf12, Y1);
99  *eqns[Y2.name()] -= dmdtf12*Y1 + fvm::Sp(dmdtf21, Y2);
100  }
101  }
102 }
103 
104 
105 template<class BasePhaseSystem>
107 (
108  const phaseSystem::dmidtfTable& dmidtfs,
110 ) const
111 {
112  forAllConstIter(phaseSystem::dmidtfTable, dmidtfs, dmidtfIter)
113  {
114  const phaseInterface interface(*this, dmidtfIter.key());
115 
116  const phaseModel& phase1 = interface.phase1();
117  const phaseModel& phase2 = interface.phase2();
118 
119  forAllConstIter(HashPtrTable<volScalarField>, *dmidtfIter(), dmidtfJter)
120  {
121  const word& member = dmidtfJter.key();
122 
123  const volScalarField& dmidtf = *dmidtfJter();
124 
125  if (!phase1.pure())
126  {
127  const volScalarField& Y1 = phase1.Y(member);
128  *eqns[Y1.name()] += dmidtf;
129  }
130 
131  if (!phase2.pure())
132  {
133  const volScalarField& Y2 = phase2.Y(member);
134  *eqns[Y2.name()] -= dmidtf;
135  }
136  }
137  }
138 }
139 
140 
141 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
142 
143 template<class BasePhaseSystem>
145 (
146  const fvMesh& mesh
147 )
148 :
149  BasePhaseSystem(mesh)
150 {
151  this->generateInterfacialModels(phaseTransferModels_);
152 
154  (
156  phaseTransferModels_,
157  phaseTransferModelIter
158  )
159  {
160  const phaseInterface& interface = phaseTransferModelIter()->interface();
161 
162  this->template validateMassTransfer<phaseTransferModel>(interface);
163 
164  if (phaseTransferModelIter()->mixture())
165  {
166  dmdtfs_.insert
167  (
168  interface,
169  new volScalarField
170  (
171  IOobject
172  (
174  (
175  "phaseTransfer:dmdtf",
176  interface.name()
177  ),
178  this->mesh().time().name(),
179  this->mesh()
180  ),
181  this->mesh(),
183  )
184  );
185 
186  d2mdtdpfs_.insert
187  (
188  interface,
189  new volScalarField
190  (
191  IOobject
192  (
194  (
195  "phaseTransfer:d2mdtdpf",
196  interface.name()
197  ),
198  this->mesh().time().name(),
199  this->mesh()
200  ),
201  this->mesh(),
203  )
204  );
205  }
206 
207  dmidtfs_.insert(interface, new HashPtrTable<volScalarField>());
208 
209  const hashedWordList species(phaseTransferModelIter()->species());
210 
211  forAllConstIter(hashedWordList, species, specieIter)
212  {
213  const word& specie = *specieIter;
214 
215  dmidtfs_[interface]->insert
216  (
217  specie,
218  new volScalarField
219  (
220  IOobject
221  (
223  (
225  (
226  "phaseTransfer:dmidtf",
227  specie
228  ),
229  interface.name()
230  ),
231  this->mesh().time().name(),
232  this->mesh()
233  ),
234  this->mesh(),
236  )
237  );
238  }
239  }
240 }
241 
242 
243 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
244 
245 template<class BasePhaseSystem>
247 {}
248 
249 
250 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
251 
252 template<class BasePhaseSystem>
255 (
256  const phaseInterfaceKey& key
257 ) const
258 {
259  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
260 
261  if (phaseTransferModels_.found(key))
262  {
263  if (phaseTransferModels_[key]->mixture())
264  {
265  tDmdtf.ref() += *dmdtfs_[key];
266  }
267 
269  (
271  *dmidtfs_[key],
272  dmidtfIter
273  )
274  {
275  tDmdtf.ref() += *dmidtfIter();
276  }
277  }
278 
279  return tDmdtf;
280 }
281 
282 
283 template<class BasePhaseSystem>
286 {
287  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
288 
289  autoPtr<phaseSystem::dmdtfTable> totalDmdtfsPtr = this->totalDmdtfs();
290  const phaseSystem::dmdtfTable& totalDmdtfs = totalDmdtfsPtr();
291 
292  forAllConstIter(phaseSystem::dmdtfTable, totalDmdtfs, totalDmdtfIter)
293  {
294  const phaseInterface interface(*this, totalDmdtfIter.key());
295 
296  addField(interface.phase1(), "dmdt", *totalDmdtfIter(), dmdts);
297  addField(interface.phase2(), "dmdt", - *totalDmdtfIter(), dmdts);
298  }
299 
300  return dmdts;
301 }
302 
303 
304 template<class BasePhaseSystem>
307 {
308  PtrList<volScalarField> d2mdtdps(BasePhaseSystem::d2mdtdps());
309 
310  forAllConstIter(phaseSystem::dmdtfTable, d2mdtdpfs_, d2mdtdpfIter)
311  {
312  const phaseInterface interface(*this, d2mdtdpfIter.key());
313 
314  addField(interface.phase1(), "d2mdtdp", *d2mdtdpfIter(), d2mdtdps);
315  addField(interface.phase2(), "d2mdtdp", - *d2mdtdpfIter(), d2mdtdps);
316  }
317 
318  return d2mdtdps;
319 }
320 
321 
322 template<class BasePhaseSystem>
325 {
327  BasePhaseSystem::momentumTransfer();
328 
329  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
330 
331  this->addDmdtUfs(totalDmdtfs(), eqns);
332 
333  return eqnsPtr;
334 }
335 
336 
337 template<class BasePhaseSystem>
340 {
342  BasePhaseSystem::momentumTransferf();
343 
344  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
345 
346  this->addDmdtUfs(totalDmdtfs(), eqns);
347 
348  return eqnsPtr;
349 }
350 
351 
352 template<class BasePhaseSystem>
355 {
357  BasePhaseSystem::heatTransfer();
358 
359  phaseSystem::heatTransferTable& eqns = eqnsPtr();
360 
361  this->addDmdtHefs(dmdtfs_, eqns);
362  this->addDmidtHefs(dmidtfs_, eqns);
363 
364  return eqnsPtr;
365 }
366 
367 
368 template<class BasePhaseSystem>
371 {
373  (
375  );
376 
377  phaseSystem::specieTransferTable& eqns = eqnsPtr();
378 
379  // Create a mass transfer matrix for each species of each phase
380  forAll(this->phaseModels_, phasei)
381  {
382  const phaseModel& phase = this->phaseModels_[phasei];
383 
384  const PtrList<volScalarField>& Yi = phase.Y();
385 
386  forAll(Yi, i)
387  {
388  eqns.insert
389  (
390  Yi[i].name(),
391  new fvScalarMatrix(Yi[i], dimMass/dimTime)
392  );
393  }
394  }
395 
396  this->addDmdtYfs(dmdtfs_, eqns);
397  this->addDmidtYf(dmidtfs_, eqns);
398 
399  return eqnsPtr;
400 }
401 
402 
403 template<class BasePhaseSystem>
405 {
407 
408  // Reset all the mass transfer rates to zero
410  (
412  phaseTransferModels_,
413  phaseTransferModelIter
414  )
415  {
416  const phaseInterface& interface = phaseTransferModelIter()->interface();
417 
418  if (phaseTransferModelIter()->mixture())
419  {
420  *dmdtfs_[interface] = Zero;
421  *d2mdtdpfs_[interface] = Zero;
422  }
423 
424  const hashedWordList species(phaseTransferModelIter()->species());
425 
426  forAllConstIter(hashedWordList, species, specieIter)
427  {
428  const word& specie = *specieIter;
429 
430  *(*dmidtfs_[interface])[specie] = Zero;
431  }
432  }
433 
434  // Evaluate the models and sum the results into the mass transfer tables
435  forAllIter
436  (
438  phaseTransferModels_,
439  phaseTransferModelIter
440  )
441  {
442  const phaseInterface& interface = phaseTransferModelIter()->interface();
443 
444  if (phaseTransferModelIter()->mixture())
445  {
446  *dmdtfs_[interface] += phaseTransferModelIter()->dmdtf();
447  *d2mdtdpfs_[interface] += phaseTransferModelIter()->d2mdtdpf();
448  }
449 
450  const HashPtrTable<volScalarField> dmidtf
451  (
452  phaseTransferModelIter()->dmidtf()
453  );
454 
455  forAllConstIter(HashPtrTable<volScalarField>, dmidtf, dmidtfIter)
456  {
457  *(*dmidtfs_[interface])[dmidtfIter.key()] += *dmidtfIter();
458  }
459  }
460 }
461 
462 
463 template<class BasePhaseSystem>
465 {
466  if (BasePhaseSystem::read())
467  {
468  bool readOK = true;
469 
470  // Models ...
471 
472  return readOK;
473  }
474  else
475  {
476  return false;
477  }
478 }
479 
480 
481 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Generic GeometricField class.
bool insert(const Key &, const T &newElmt)
Insert a new hashedEntry.
Definition: HashTableI.H:80
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:149
const word & name() const
Return name.
Definition: IOobject.H:310
static word groupName(Name name, const word &group)
Class which models non-thermally-coupled or weakly thermally coupled mass transfers.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
virtual void correct()
Correct the mass transfer rates.
void addDmidtYf(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from specie mass transfers.
virtual tmp< volScalarField > dmdtf(const phaseInterfaceKey &key) const
Return the mass transfer rate for an interface.
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
PhaseTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
virtual ~PhaseTransferPhaseSystem()
Destructor.
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
void addDmdtYfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from bulk mass transfers.
virtual PtrList< volScalarField > d2mdtdps() const
Return the mass transfer pressure linearisation coefficients for.
virtual bool read()
Read base phaseProperties dictionary.
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
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: autoPtr.H:51
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:101
A wordList with hashed indices for faster lookup by name.
Word-pair based class used for keying interface models in hash tables.
Class to represent an interface between phases. Derivations can further specify the configuration of ...
const phaseModel & phase1() const
Return phase 1.
const phaseModel & phase2() const
Return phase 2.
virtual bool pure() const =0
Return whether the phase is pure (i.e., not multi-component)
virtual const PtrList< volScalarField > & Y() const =0
Return the species mass fractions.
Base class of the thermophysical property types.
Definition: specie.H:66
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
Calculate the matrix for implicit and explicit sources.
void correct(const RdeltaTType &rDeltaT, const RhoType &rho, volScalarField &psi, const surfaceScalarField &phiCorr, const SpType &Sp, const SuType &Su)
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
tmp< VolField< Type > > Sp(const volScalarField &sp, const VolField< Type > &vf)
Definition: fvcSup.C:67
void addField(const Group &group, const word &name, tmp< GeoField > field, PtrList< GeoField > &fieldList)
static const zero Zero
Definition: zero.H:97
const dimensionSet dimPressure
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
const dimensionSet dimTime
dimensionedScalar negPart(const dimensionedScalar &ds)
const dimensionSet dimDensity
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar posPart(const dimensionedScalar &ds)