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-2022 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,
79  phaseSystem::specieTransferTable& eqns
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,
109  phaseSystem::specieTransferTable& eqns
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  (
155  phaseTransferModelTable,
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().timeName(),
179  this->mesh()
180  ),
181  this->mesh(),
183  )
184  );
185  }
186 
187  dmidtfs_.insert(interface, new HashPtrTable<volScalarField>());
188 
189  const hashedWordList species(phaseTransferModelIter()->species());
190 
191  forAllConstIter(hashedWordList, species, specieIter)
192  {
193  const word& specie = *specieIter;
194 
195  dmidtfs_[interface]->insert
196  (
197  specie,
198  new volScalarField
199  (
200  IOobject
201  (
203  (
205  (
206  "phaseTransfer:dmidtf",
207  specie
208  ),
209  interface.name()
210  ),
211  this->mesh().time().timeName(),
212  this->mesh()
213  ),
214  this->mesh(),
216  )
217  );
218  }
219  }
220 }
221 
222 
223 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
224 
225 template<class BasePhaseSystem>
227 {}
228 
229 
230 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
231 
232 template<class BasePhaseSystem>
235 (
236  const phaseInterfaceKey& key
237 ) const
238 {
239  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
240 
241  if (phaseTransferModels_.found(key))
242  {
243  if (phaseTransferModels_[key]->mixture())
244  {
245  tDmdtf.ref() += *dmdtfs_[key];
246  }
247 
249  (
250  HashPtrTable<volScalarField>,
251  *dmidtfs_[key],
252  dmidtfIter
253  )
254  {
255  tDmdtf.ref() += *dmidtfIter();
256  }
257  }
258 
259  return tDmdtf;
260 }
261 
262 
263 template<class BasePhaseSystem>
266 {
267  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
268 
269  autoPtr<phaseSystem::dmdtfTable> totalDmdtfsPtr = this->totalDmdtfs();
270  const phaseSystem::dmdtfTable& totalDmdtfs = totalDmdtfsPtr();
271 
272  forAllConstIter(phaseSystem::dmdtfTable, totalDmdtfs, totalDmdtfIter)
273  {
274  const phaseInterface interface(*this, totalDmdtfIter.key());
275 
276  addField(interface.phase1(), "dmdt", *totalDmdtfIter(), dmdts);
277  addField(interface.phase2(), "dmdt", - *totalDmdtfIter(), dmdts);
278  }
279 
280  return dmdts;
281 }
282 
283 
284 template<class BasePhaseSystem>
287 {
288  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
290 
291  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
292 
293  this->addDmdtUfs(totalDmdtfs(), eqns);
294 
295  return eqnsPtr;
296 }
297 
298 
299 template<class BasePhaseSystem>
302 {
303  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
304  BasePhaseSystem::momentumTransferf();
305 
306  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
307 
308  this->addDmdtUfs(totalDmdtfs(), eqns);
309 
310  return eqnsPtr;
311 }
312 
313 
314 template<class BasePhaseSystem>
317 {
318  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
319  BasePhaseSystem::heatTransfer();
320 
321  phaseSystem::heatTransferTable& eqns = eqnsPtr();
322 
323  this->addDmdtHefs(dmdtfs_, eqns);
324  this->addDmidtHefs(dmidtfs_, eqns);
325 
326  return eqnsPtr;
327 }
328 
329 
330 template<class BasePhaseSystem>
333 {
334  autoPtr<phaseSystem::specieTransferTable> eqnsPtr
335  (
337  );
338 
339  phaseSystem::specieTransferTable& eqns = eqnsPtr();
340 
341  // Create a mass transfer matrix for each species of each phase
342  forAll(this->phaseModels_, phasei)
343  {
344  const phaseModel& phase = this->phaseModels_[phasei];
345 
346  const PtrList<volScalarField>& Yi = phase.Y();
347 
348  forAll(Yi, i)
349  {
350  eqns.insert
351  (
352  Yi[i].name(),
353  new fvScalarMatrix(Yi[i], dimMass/dimTime)
354  );
355  }
356  }
357 
358  this->addDmdtYfs(dmdtfs_, eqns);
359  this->addDmidtYf(dmidtfs_, eqns);
360 
361  return eqnsPtr;
362 }
363 
364 
365 template<class BasePhaseSystem>
367 {
369 
370  // Reset all the mass transfer rates to zero
372  (
373  phaseTransferModelTable,
374  phaseTransferModels_,
375  phaseTransferModelIter
376  )
377  {
378  const phaseInterface& interface = phaseTransferModelIter()->interface();
379 
380  if (phaseTransferModelIter()->mixture())
381  {
382  *dmdtfs_[interface] = Zero;
383  }
384 
385  const hashedWordList species(phaseTransferModelIter()->species());
386 
387  forAllConstIter(hashedWordList, species, specieIter)
388  {
389  const word& specie = *specieIter;
390 
391  *(*dmidtfs_[interface])[specie] = Zero;
392  }
393  }
394 
395  // Evaluate the models and sum the results into the mass transfer tables
396  forAllIter
397  (
398  phaseTransferModelTable,
399  phaseTransferModels_,
400  phaseTransferModelIter
401  )
402  {
403  const phaseInterface& interface = phaseTransferModelIter()->interface();
404 
405  if (phaseTransferModelIter()->mixture())
406  {
407  *dmdtfs_[interface] += phaseTransferModelIter()->dmdtf();
408  }
409 
410  const HashPtrTable<volScalarField> dmidtf
411  (
412  phaseTransferModelIter()->dmidtf()
413  );
414 
415  forAllConstIter(HashPtrTable<volScalarField>, dmidtf, dmidtfIter)
416  {
417  *(*dmidtfs_[interface])[dmidtfIter.key()] += *dmidtfIter();
418  }
419  }
420 }
421 
422 
423 template<class BasePhaseSystem>
425 {
426  if (BasePhaseSystem::read())
427  {
428  bool readOK = true;
429 
430  // Models ...
431 
432  return readOK;
433  }
434  else
435  {
436  return false;
437  }
438 }
439 
440 
441 // ************************************************************************* //
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
PhaseTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
label phasei
Definition: pEqn.H:27
#define forAllConstIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:477
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-10/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:76
virtual autoPtr< phaseSystem::heatTransferTable > heatTransfer() const
Return the heat transfer matrices.
#define forAllIter(Container, container, iter)
Iterate across all elements in the container object of type.
Definition: UList.H:459
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransferf()
Return the momentum transfer matrices for the face-based algorithm.
Class which models non-thermally-coupled or weakly thermally coupled mass transfers.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:78
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:80
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
fvMesh & mesh
dimensionedScalar posPart(const dimensionedScalar &ds)
virtual void correct()
Correct the mass transfer rates.
void addDmdtHefs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from bulk mass transfers.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:58
void addDmidtHefs(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie mass transfers.
const dimensionSet dimTime
virtual bool read()
Read base phaseProperties dictionary.
static word groupName(Name name, const word &group)
const dimensionSet dimDensity
static const zero Zero
Definition: zero.H:97
void addDmdtYfs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::specieTransferTable &eqns) const
Add specie transfer terms which result from bulk mass transfers.
HashPtrTable< volScalarField, phaseInterfaceKey, phaseInterfaceKey::hash > dmdtfTable
Definition: phaseSystem.H:93
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
const dimensionSet dimMass
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
phaseSystem::momentumTransferTable & momentumTransfer(momentumTransferPtr())
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
virtual autoPtr< phaseSystem::specieTransferTable > specieTransfer() const
Return the specie transfer matrices.
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.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:70
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
A class for managing temporary objects.
Definition: PtrList.H:53
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)
virtual ~PhaseTransferPhaseSystem()
Destructor.