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-2020 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 phasePair& pair =
50  this->phasePairs_[phaseTransferModelIter.key()];
51 
52  totalDmdtfs.insert(pair, phaseSystem::dmdtf(pair).ptr());
53 
54  if (phaseTransferModelIter()->mixture())
55  {
56  *totalDmdtfs[pair] += *dmdtfs_[pair];
57  }
58 
60  (
61  HashPtrTable<volScalarField>,
62  *dmidtfs_[pair],
63  dmidtfIter
64  )
65  {
66  *totalDmdtfs[pair] += *dmidtfIter();
67  }
68  }
69 
70  return totalDmdtfsPtr;
71 }
72 
73 
74 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
75 
76 template<class BasePhaseSystem>
78 (
79  const phaseSystem::dmdtfTable& dmdtfs,
80  phaseSystem::specieTransferTable& eqns
81 ) const
82 {
83  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
84  {
85  const phasePairKey& key = dmdtfIter.key();
86  const phasePair& pair(this->phasePairs_[key]);
87 
88  const volScalarField dmdtf(Pair<word>::compare(pair, key)**dmdtfIter());
89  const volScalarField dmdtf12(negPart(dmdtf));
90  const volScalarField dmdtf21(posPart(dmdtf));
91 
92  const phaseModel& phase1 = pair.phase1();
93  const phaseModel& phase2 = pair.phase2();
94 
95  forAll(phase1.Y(), Yi1)
96  {
97  const volScalarField& Y1 = phase1.Y()[Yi1];
98  const volScalarField& Y2 = phase2.Y(Y1.member());
99 
100  *eqns[Y1.name()] += dmdtf21*Y2 + fvm::Sp(dmdtf12, Y1);
101  *eqns[Y2.name()] -= dmdtf12*Y1 + fvm::Sp(dmdtf21, Y2);
102  }
103  }
104 }
105 
106 
107 template<class BasePhaseSystem>
109 (
110  const phaseSystem::dmidtfTable& dmidtfs,
111  phaseSystem::specieTransferTable& eqns
112 ) const
113 {
114  forAllConstIter(phaseSystem::dmidtfTable, dmidtfs, dmidtfIter)
115  {
116  const phasePairKey& key = dmidtfIter.key();
117  const phasePair& pair(this->phasePairs_[key]);
118 
119  const phaseModel& phase1 = pair.phase1();
120  const phaseModel& phase2 = pair.phase2();
121 
122  forAllConstIter(HashPtrTable<volScalarField>, *dmidtfIter(), dmidtfJter)
123  {
124  const word& member = dmidtfJter.key();
125 
126  const volScalarField dmidtf
127  (
128  Pair<word>::compare(pair, key)**dmidtfJter()
129  );
130 
131  if (!phase1.pure())
132  {
133  const volScalarField& Y1 = phase1.Y(member);
134  *eqns[Y1.name()] += dmidtf;
135  }
136 
137  if (!phase2.pure())
138  {
139  const volScalarField& Y2 = phase2.Y(member);
140  *eqns[Y2.name()] -= dmidtf;
141  }
142  }
143  }
144 }
145 
146 
147 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
148 
149 template<class BasePhaseSystem>
151 (
152  const fvMesh& mesh
153 )
154 :
155  BasePhaseSystem(mesh)
156 {
157  this->generatePairsAndSubModels
158  (
159  "phaseTransfer",
160  phaseTransferModels_,
161  false
162  );
163 
165  (
166  phaseTransferModelTable,
167  phaseTransferModels_,
168  phaseTransferModelIter
169  )
170  {
171  const phasePair& pair = this->phasePairs_[phaseTransferModelIter.key()];
172 
173  if (phaseTransferModelIter()->mixture())
174  {
175  dmdtfs_.insert
176  (
177  pair,
178  new volScalarField
179  (
180  IOobject
181  (
183  (
184  "phaseTransfer:dmdtf",
185  pair.name()
186  ),
187  this->mesh().time().timeName(),
188  this->mesh()
189  ),
190  this->mesh(),
192  )
193  );
194  }
195 
196  dmidtfs_.insert(pair, new HashPtrTable<volScalarField>());
197 
198  const hashedWordList species(phaseTransferModelIter()->species());
199 
200  forAllConstIter(hashedWordList, species, specieIter)
201  {
202  const word& specie = *specieIter;
203 
204  dmidtfs_[pair]->insert
205  (
206  specie,
207  new volScalarField
208  (
209  IOobject
210  (
212  (
214  (
215  "phaseTransfer:dmidtf",
216  specie
217  ),
218  pair.name()
219  ),
220  this->mesh().time().timeName(),
221  this->mesh()
222  ),
223  this->mesh(),
225  )
226  );
227  }
228  }
229 }
230 
231 
232 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
233 
234 template<class BasePhaseSystem>
236 {}
237 
238 
239 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
240 
241 template<class BasePhaseSystem>
244 (
245  const phasePairKey& key
246 ) const
247 {
248  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
249 
250  if (phaseTransferModels_.found(key))
251  {
252  const label dmdtfSign(Pair<word>::compare(this->phasePairs_[key], key));
253 
254  if (phaseTransferModels_[key]->mixture())
255  {
256  tDmdtf.ref() += dmdtfSign**dmdtfs_[key];
257  }
258 
260  (
261  HashPtrTable<volScalarField>,
262  *dmidtfs_[key],
263  dmidtfIter
264  )
265  {
266  tDmdtf.ref() += dmdtfSign**dmidtfIter();
267  }
268  }
269 
270  return tDmdtf;
271 }
272 
273 
274 template<class BasePhaseSystem>
277 {
278  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
279 
280  autoPtr<phaseSystem::dmdtfTable> totalDmdtfsPtr = this->totalDmdtfs();
281  const phaseSystem::dmdtfTable& totalDmdtfs = totalDmdtfsPtr();
282 
283  forAllConstIter(phaseSystem::dmdtfTable, totalDmdtfs, totalDmdtfIter)
284  {
285  const phasePair& pair = this->phasePairs_[totalDmdtfIter.key()];
286 
287  addField(pair.phase1(), "dmdt", *totalDmdtfIter(), dmdts);
288  addField(pair.phase2(), "dmdt", - *totalDmdtfIter(), dmdts);
289  }
290 
291  return dmdts;
292 }
293 
294 
295 template<class BasePhaseSystem>
298 {
299  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
301 
302  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
303 
304  this->addDmdtUfs(totalDmdtfs(), eqns);
305 
306  return eqnsPtr;
307 }
308 
309 
310 template<class BasePhaseSystem>
313 {
314  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
315  BasePhaseSystem::momentumTransferf();
316 
317  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
318 
319  this->addDmdtUfs(totalDmdtfs(), eqns);
320 
321  return eqnsPtr;
322 }
323 
324 
325 template<class BasePhaseSystem>
328 {
329  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
330  BasePhaseSystem::heatTransfer();
331 
332  phaseSystem::heatTransferTable& eqns = eqnsPtr();
333 
334  this->addDmdtHefs(dmdtfs_, eqns);
335  this->addDmidtHef(dmidtfs_, eqns);
336 
337  return eqnsPtr;
338 }
339 
340 
341 template<class BasePhaseSystem>
344 {
345  autoPtr<phaseSystem::specieTransferTable> eqnsPtr
346  (
348  );
349 
350  phaseSystem::specieTransferTable& eqns = eqnsPtr();
351 
352  // Create a mass transfer matrix for each species of each phase
353  forAll(this->phaseModels_, phasei)
354  {
355  const phaseModel& phase = this->phaseModels_[phasei];
356 
357  const PtrList<volScalarField>& Yi = phase.Y();
358 
359  forAll(Yi, i)
360  {
361  eqns.insert
362  (
363  Yi[i].name(),
364  new fvScalarMatrix(Yi[i], dimMass/dimTime)
365  );
366  }
367  }
368 
369  this->addDmdtYfs(dmdtfs_, eqns);
370  this->addDmidtYf(dmidtfs_, eqns);
371 
372  return eqnsPtr;
373 }
374 
375 
376 template<class BasePhaseSystem>
378 {
380 
381  // Reset all the mass transfer rates to zero
383  (
384  phaseTransferModelTable,
385  phaseTransferModels_,
386  phaseTransferModelIter
387  )
388  {
389  const phasePair& pair = this->phasePairs_[phaseTransferModelIter.key()];
390 
391  if (phaseTransferModelIter()->mixture())
392  {
393  *dmdtfs_[pair] = Zero;
394  }
395 
396  const hashedWordList species(phaseTransferModelIter()->species());
397 
398  forAllConstIter(hashedWordList, species, specieIter)
399  {
400  const word& specie = *specieIter;
401 
402  *(*dmidtfs_[pair])[specie] = Zero;
403  }
404  }
405 
406  // Evaluate the models and sum the results into the mass transfer tables
407  forAllIter
408  (
409  phaseTransferModelTable,
410  phaseTransferModels_,
411  phaseTransferModelIter
412  )
413  {
414  const phasePair& pair = this->phasePairs_[phaseTransferModelIter.key()];
415 
416  if (phaseTransferModelIter()->mixture())
417  {
418  *dmdtfs_[pair] += phaseTransferModelIter()->dmdtf();
419  }
420 
421  const HashPtrTable<volScalarField> dmidtf
422  (
423  phaseTransferModelIter()->dmidtf()
424  );
425 
426  forAllConstIter(HashPtrTable<volScalarField>, dmidtf, dmidtfIter)
427  {
428  *(*dmidtfs_[pair])[dmidtfIter.key()] += *dmidtfIter();
429  }
430  }
431 }
432 
433 
434 template<class BasePhaseSystem>
436 {
437  if (BasePhaseSystem::read())
438  {
439  bool readOK = true;
440 
441  // Models ...
442 
443  return readOK;
444  }
445  else
446  {
447  return false;
448  }
449 }
450 
451 
452 // ************************************************************************* //
void addDmidtHef(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie mass transfers.
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
PhaseTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
label phasei
Definition: pEqn.H:27
HashPtrTable< fvVectorMatrix > momentumTransferTable
Definition: phaseSystem.H:75
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.
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:142
Class which models non-thermally-coupled or weakly thermally coupled mass transfers.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:77
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:79
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
dimensionedScalar posPart(const dimensionedScalar &ds)
void addDmdtHefs(const phaseSystem::dmdtfTable &dmdtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from bulk mass transfers.
virtual void correct()
Correct the mass transfer rates.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:57
dynamicFvMesh & mesh
virtual bool read()
Read base phaseProperties dictionary.
static word groupName(Name name, const word &group)
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtfTable
Definition: phaseSystem.H:91
virtual tmp< volScalarField > dmdtf(const phasePairKey &key) const
Return the mass transfer rate for an interface.
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.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-8/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
virtual PtrList< volScalarField > dmdts() const
Return the mass transfer rates for each phase.
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
phaseSystem::momentumTransferTable & momentumTransfer(momentumTransferPtr())
const dimensionSet dimDensity
virtual autoPtr< phaseSystem::momentumTransferTable > momentumTransfer()
Return the momentum transfer matrices for the cell-based algorithm.
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.
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
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: PtrList.H:53
phaseChangeTwoPhaseMixture & mixture
Definition: createFields.H:38
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)
virtual ~PhaseTransferPhaseSystem()
Destructor.
zeroField Sp
Definition: alphaSuSp.H:2