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-2021 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  this->template validateMassTransfer<phaseTransferModel>(pair);
174 
175  if (phaseTransferModelIter()->mixture())
176  {
177  dmdtfs_.insert
178  (
179  pair,
180  new volScalarField
181  (
182  IOobject
183  (
185  (
186  "phaseTransfer:dmdtf",
187  pair.name()
188  ),
189  this->mesh().time().timeName(),
190  this->mesh()
191  ),
192  this->mesh(),
194  )
195  );
196  }
197 
198  dmidtfs_.insert(pair, new HashPtrTable<volScalarField>());
199 
200  const hashedWordList species(phaseTransferModelIter()->species());
201 
202  forAllConstIter(hashedWordList, species, specieIter)
203  {
204  const word& specie = *specieIter;
205 
206  dmidtfs_[pair]->insert
207  (
208  specie,
209  new volScalarField
210  (
211  IOobject
212  (
214  (
216  (
217  "phaseTransfer:dmidtf",
218  specie
219  ),
220  pair.name()
221  ),
222  this->mesh().time().timeName(),
223  this->mesh()
224  ),
225  this->mesh(),
227  )
228  );
229  }
230  }
231 }
232 
233 
234 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
235 
236 template<class BasePhaseSystem>
238 {}
239 
240 
241 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
242 
243 template<class BasePhaseSystem>
246 (
247  const phasePairKey& key
248 ) const
249 {
250  tmp<volScalarField> tDmdtf = BasePhaseSystem::dmdtf(key);
251 
252  if (phaseTransferModels_.found(key))
253  {
254  const label dmdtfSign(Pair<word>::compare(this->phasePairs_[key], key));
255 
256  if (phaseTransferModels_[key]->mixture())
257  {
258  tDmdtf.ref() += dmdtfSign**dmdtfs_[key];
259  }
260 
262  (
263  HashPtrTable<volScalarField>,
264  *dmidtfs_[key],
265  dmidtfIter
266  )
267  {
268  tDmdtf.ref() += dmdtfSign**dmidtfIter();
269  }
270  }
271 
272  return tDmdtf;
273 }
274 
275 
276 template<class BasePhaseSystem>
279 {
280  PtrList<volScalarField> dmdts(BasePhaseSystem::dmdts());
281 
282  autoPtr<phaseSystem::dmdtfTable> totalDmdtfsPtr = this->totalDmdtfs();
283  const phaseSystem::dmdtfTable& totalDmdtfs = totalDmdtfsPtr();
284 
285  forAllConstIter(phaseSystem::dmdtfTable, totalDmdtfs, totalDmdtfIter)
286  {
287  const phasePair& pair = this->phasePairs_[totalDmdtfIter.key()];
288 
289  addField(pair.phase1(), "dmdt", *totalDmdtfIter(), dmdts);
290  addField(pair.phase2(), "dmdt", - *totalDmdtfIter(), dmdts);
291  }
292 
293  return dmdts;
294 }
295 
296 
297 template<class BasePhaseSystem>
300 {
301  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
303 
304  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
305 
306  this->addDmdtUfs(totalDmdtfs(), eqns);
307 
308  return eqnsPtr;
309 }
310 
311 
312 template<class BasePhaseSystem>
315 {
316  autoPtr<phaseSystem::momentumTransferTable> eqnsPtr =
317  BasePhaseSystem::momentumTransferf();
318 
319  phaseSystem::momentumTransferTable& eqns = eqnsPtr();
320 
321  this->addDmdtUfs(totalDmdtfs(), eqns);
322 
323  return eqnsPtr;
324 }
325 
326 
327 template<class BasePhaseSystem>
330 {
331  autoPtr<phaseSystem::heatTransferTable> eqnsPtr =
332  BasePhaseSystem::heatTransfer();
333 
334  phaseSystem::heatTransferTable& eqns = eqnsPtr();
335 
336  this->addDmdtHefs(dmdtfs_, eqns);
337  this->addDmidtHefs(dmidtfs_, eqns);
338 
339  return eqnsPtr;
340 }
341 
342 
343 template<class BasePhaseSystem>
346 {
347  autoPtr<phaseSystem::specieTransferTable> eqnsPtr
348  (
350  );
351 
352  phaseSystem::specieTransferTable& eqns = eqnsPtr();
353 
354  // Create a mass transfer matrix for each species of each phase
355  forAll(this->phaseModels_, phasei)
356  {
357  const phaseModel& phase = this->phaseModels_[phasei];
358 
359  const PtrList<volScalarField>& Yi = phase.Y();
360 
361  forAll(Yi, i)
362  {
363  eqns.insert
364  (
365  Yi[i].name(),
366  new fvScalarMatrix(Yi[i], dimMass/dimTime)
367  );
368  }
369  }
370 
371  this->addDmdtYfs(dmdtfs_, eqns);
372  this->addDmidtYf(dmidtfs_, eqns);
373 
374  return eqnsPtr;
375 }
376 
377 
378 template<class BasePhaseSystem>
380 {
382 
383  // Reset all the mass transfer rates to zero
385  (
386  phaseTransferModelTable,
387  phaseTransferModels_,
388  phaseTransferModelIter
389  )
390  {
391  const phasePair& pair = this->phasePairs_[phaseTransferModelIter.key()];
392 
393  if (phaseTransferModelIter()->mixture())
394  {
395  *dmdtfs_[pair] = Zero;
396  }
397 
398  const hashedWordList species(phaseTransferModelIter()->species());
399 
400  forAllConstIter(hashedWordList, species, specieIter)
401  {
402  const word& specie = *specieIter;
403 
404  *(*dmidtfs_[pair])[specie] = Zero;
405  }
406  }
407 
408  // Evaluate the models and sum the results into the mass transfer tables
409  forAllIter
410  (
411  phaseTransferModelTable,
412  phaseTransferModels_,
413  phaseTransferModelIter
414  )
415  {
416  const phasePair& pair = this->phasePairs_[phaseTransferModelIter.key()];
417 
418  if (phaseTransferModelIter()->mixture())
419  {
420  *dmdtfs_[pair] += phaseTransferModelIter()->dmdtf();
421  }
422 
423  const HashPtrTable<volScalarField> dmidtf
424  (
425  phaseTransferModelIter()->dmidtf()
426  );
427 
428  forAllConstIter(HashPtrTable<volScalarField>, dmidtf, dmidtfIter)
429  {
430  *(*dmidtfs_[pair])[dmidtfIter.key()] += *dmidtfIter();
431  }
432  }
433 }
434 
435 
436 template<class BasePhaseSystem>
438 {
439  if (BasePhaseSystem::read())
440  {
441  bool readOK = true;
442 
443  // Models ...
444 
445  return readOK;
446  }
447  else
448  {
449  return false;
450  }
451 }
452 
453 
454 // ************************************************************************* //
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:78
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:153
Class which models non-thermally-coupled or weakly thermally coupled mass transfers.
HashPtrTable< fvScalarMatrix > heatTransferTable
Definition: phaseSystem.H:80
HashPtrTable< fvScalarMatrix > specieTransferTable
Definition: phaseSystem.H:82
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
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:57
void addDmidtHefs(const phaseSystem::dmidtfTable &dmidtfs, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie mass transfers.
const dimensionSet dimTime
dynamicFvMesh & mesh
virtual bool read()
Read base phaseProperties dictionary.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-9/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
static word groupName(Name name, const word &group)
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtfTable
Definition: phaseSystem.H:94
virtual tmp< volScalarField > dmdtf(const phasePairKey &key) const
Return the mass transfer rate for an interface.
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.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
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.
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.