HeatTransferPhaseSystem.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) 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 "fvmSup.H"
28 #include "rhoReactionThermo.H"
29 
30 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
31 
32 template<class BasePhaseSystem>
34 (
35  const phaseSystem::dmdtfTable& dmdtfs,
36  phaseSystem::heatTransferTable& eqns
37 ) const
38 {
39  // Loop the pairs
40  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
41  {
42  const phasePairKey& key = dmdtfIter.key();
43  const phasePair& pair(this->phasePairs_[key]);
44 
45  const volScalarField dmdtf(Pair<word>::compare(pair, key)**dmdtfIter());
46  const volScalarField dmdtf21(posPart(dmdtf));
47  const volScalarField dmdtf12(negPart(dmdtf));
48 
49  const phaseModel& phase1 = pair.phase1();
50  const phaseModel& phase2 = pair.phase2();
51  const rhoThermo& thermo1 = phase1.thermo();
52  const rhoThermo& thermo2 = phase2.thermo();
53  const volScalarField& he1 = thermo1.he();
54  const volScalarField& he2 = thermo2.he();
55  const volScalarField hs1(thermo1.hs());
56  const volScalarField hs2(thermo2.hs());
57  const volScalarField K1(phase1.K());
58  const volScalarField K2(phase2.K());
59 
60  // Transfer of sensible enthalpy within the phases
61  *eqns[phase1.name()] +=
62  dmdtf*hs1 + fvm::Sp(dmdtf12, he1) - dmdtf12*he1;
63  *eqns[phase2.name()] -=
64  dmdtf*hs2 + fvm::Sp(dmdtf21, he2) - dmdtf21*he2;
65 
66  // Transfer of sensible enthalpy between the phases
67  *eqns[phase1.name()] += dmdtf21*(hs2 - hs1);
68  *eqns[phase2.name()] -= dmdtf12*(hs1 - hs2);
69 
70  // Transfer of kinetic energy
71  *eqns[phase1.name()] += dmdtf21*K2 + dmdtf12*K1;
72  *eqns[phase2.name()] -= dmdtf12*K1 + dmdtf21*K2;
73  }
74 }
75 
76 
77 template<class BasePhaseSystem>
79 (
80  const phaseSystem::dmidtfTable& dmidtfs,
81  phaseSystem::heatTransferTable& eqns
82 ) const
83 {
84  static const dimensionedScalar one(dimless, 1);
85 
86  // Loop the pairs
87  forAllConstIter(phaseSystem::dmidtfTable, dmidtfs, dmidtfIter)
88  {
89  const phasePairKey& key = dmidtfIter.key();
90  const phasePair& pair(this->phasePairs_[key]);
91 
92  const phaseModel& phase1 = pair.phase1();
93  const phaseModel& phase2 = pair.phase2();
94  const rhoThermo& thermo1 = phase1.thermo();
95  const rhoThermo& thermo2 = phase2.thermo();
96  const basicSpecieMixture* compositionPtr1 =
97  isA<rhoReactionThermo>(thermo1)
98  ? &refCast<const rhoReactionThermo>(thermo1).composition()
99  : static_cast<const basicSpecieMixture*>(nullptr);
100  const basicSpecieMixture* compositionPtr2 =
101  isA<rhoReactionThermo>(thermo2)
102  ? &refCast<const rhoReactionThermo>(thermo2).composition()
103  : static_cast<const basicSpecieMixture*>(nullptr);
104  const volScalarField& he1 = thermo1.he();
105  const volScalarField& he2 = thermo2.he();
106  const volScalarField hs1(thermo1.hs());
107  const volScalarField hs2(thermo2.hs());
108  const volScalarField K1(phase1.K());
109  const volScalarField K2(phase2.K());
110 
111  // Loop the species
112  forAllConstIter(HashPtrTable<volScalarField>, *dmidtfIter(), dmidtfJter)
113  {
114  const word& specie = dmidtfJter.key();
115 
116  // Mass transfer rates
117  const volScalarField dmidtf
118  (
119  Pair<word>::compare(pair, key)**dmidtfJter()
120  );
121  const volScalarField dmidtf21(posPart(dmidtf));
122  const volScalarField dmidtf12(negPart(dmidtf));
123 
124  // Specie indices
125  const label speciei1 =
126  compositionPtr1 ? compositionPtr1->species()[specie] : -1;
127  const label speciei2 =
128  compositionPtr2 ? compositionPtr2->species()[specie] : -1;
129 
130  // Enthalpies
131  const volScalarField hsi1
132  (
133  compositionPtr1
134  ? compositionPtr1->Hs(speciei1, thermo1.p(), thermo1.T())
135  : tmp<volScalarField>(hs1)
136  );
137  const volScalarField hsi2
138  (
139  compositionPtr2
140  ? compositionPtr2->Hs(speciei2, thermo2.p(), thermo2.T())
141  : tmp<volScalarField>(hs2)
142  );
143 
144  // Limited mass fractions
145  tmp<volScalarField> tYi1, tYi2;
146  if (residualY_ > 0)
147  {
148  tYi1 =
149  compositionPtr1
150  ? max(compositionPtr1->Y(speciei1), residualY_)
151  : volScalarField::New("Yi1", this->mesh(), one);
152  tYi2 =
153  compositionPtr2
154  ? max(compositionPtr2->Y(speciei2), residualY_)
155  : volScalarField::New("Yi2", this->mesh(), one);
156  }
157 
158  // Transfer of sensible enthalpy within the phases
159  *eqns[phase1.name()] += dmidtf*hsi1;
160  *eqns[phase2.name()] -= dmidtf*hsi2;
161  if (residualY_ > 0)
162  {
163  *eqns[phase1.name()] +=
164  fvm::Sp(dmidtf12/tYi1(), he1) - dmidtf12/tYi1()*he1;
165  *eqns[phase2.name()] -=
166  fvm::Sp(dmidtf21/tYi2(), he2) - dmidtf21/tYi2()*he2;
167  }
168 
169  // Transfer of sensible enthalpy between the phases
170  *eqns[phase1.name()] += dmidtf21*(hsi2 - hsi1);
171  *eqns[phase2.name()] -= dmidtf12*(hsi1 - hsi2);
172 
173  // Transfer of kinetic energy
174  *eqns[phase1.name()] += dmidtf21*K2 + dmidtf12*K1;
175  *eqns[phase2.name()] -= dmidtf12*K1 + dmidtf21*K2;
176  }
177  }
178 }
179 
180 
181 template<class BasePhaseSystem>
183 (
184  const phaseSystem::dmdtfTable& dmdtfs,
185  const phaseSystem::dmdtfTable& Tfs,
186  const latentHeatScheme scheme,
187  phaseSystem::heatTransferTable& eqns
188 ) const
189 {
190  // Loop the pairs
191  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
192  {
193  const phasePairKey& key = dmdtfIter.key();
194  const phasePair& pair(this->phasePairs_[key]);
195 
196  const volScalarField dmdtf(Pair<word>::compare(pair, key)**dmdtfIter());
197  const volScalarField dmdtf21(posPart(dmdtf));
198  const volScalarField dmdtf12(negPart(dmdtf));
199 
200  const volScalarField& Tf = *Tfs[key];
201 
202  const phaseModel& phase1 = pair.phase1();
203  const phaseModel& phase2 = pair.phase2();
204  const rhoThermo& thermo1 = phase1.thermo();
205  const rhoThermo& thermo2 = phase2.thermo();
206  const volScalarField& he1 = thermo1.he();
207  const volScalarField& he2 = thermo2.he();
208  const volScalarField K1(phase1.K());
209  const volScalarField K2(phase2.K());
210 
211  // Interface enthalpies
212  const volScalarField hsf1(thermo1.hs(thermo1.p(), Tf));
213  const volScalarField hsf2(thermo2.hs(thermo1.p(), Tf));
214 
215  // Transfer of energy from the interface into the bulk
216  switch (scheme)
217  {
218  case latentHeatScheme::symmetric:
219  {
220  *eqns[phase1.name()] += dmdtf*hsf1;
221  *eqns[phase2.name()] -= dmdtf*hsf2;
222 
223  break;
224  }
225  case latentHeatScheme::upwind:
226  {
227  // Bulk enthalpies
228  const volScalarField hs1(thermo1.hs());
229  const volScalarField hs2(thermo2.hs());
230 
231  *eqns[phase1.name()] += dmdtf21*hsf1 + dmdtf12*hs1;
232  *eqns[phase2.name()] -= dmdtf12*hsf2 + dmdtf21*hs2;
233 
234  break;
235  }
236  }
237  *eqns[phase1.name()] += fvm::Sp(dmdtf12, he1) - dmdtf12*he1;
238  *eqns[phase2.name()] -= fvm::Sp(dmdtf21, he2) - dmdtf21*he2;
239 
240  // Transfer of kinetic energy
241  *eqns[phase1.name()] += dmdtf21*K2 + dmdtf12*K1;
242  *eqns[phase2.name()] -= dmdtf12*K1 + dmdtf21*K2;
243  }
244 }
245 
246 
247 template<class BasePhaseSystem>
249 (
250  const phaseSystem::dmdtfTable& dmdtfs,
251  const phaseSystem::dmdtfTable& Tfs,
252  const scalar weight,
253  const latentHeatScheme scheme,
254  phaseSystem::heatTransferTable& eqns
255 ) const
256 {
257  // Loop the pairs
258  forAllConstIter(phaseSystem::dmdtfTable, dmdtfs, dmdtfIter)
259  {
260  const phasePairKey& key = dmdtfIter.key();
261  const phasePair& pair(this->phasePairs_[key]);
262 
263  const volScalarField dmdtf(Pair<word>::compare(pair, key)**dmdtfIter());
264  const volScalarField dmdtf21(posPart(dmdtf));
265  const volScalarField dmdtf12(negPart(dmdtf));
266 
267  const volScalarField& Tf = *Tfs[key];
268 
269  const phaseModel& phase1 = pair.phase1();
270  const phaseModel& phase2 = pair.phase2();
271 
272  // Latent heat contribution
273  const volScalarField L(this->L(pair, dmdtf, Tf, scheme));
274  *eqns[phase1.name()] += ((1 - weight)*dmdtf12 + weight*dmdtf21)*L;
275  *eqns[phase2.name()] += ((1 - weight)*dmdtf21 + weight*dmdtf12)*L;
276  }
277 }
278 
279 
280 template<class BasePhaseSystem>
282 (
283  const phaseSystem::dmdtfTable& dmdtfs,
284  const phaseSystem::dmdtfTable& Tfs,
285  const scalar weight,
286  const latentHeatScheme scheme,
287  phaseSystem::heatTransferTable& eqns
288 ) const
289 {
290  addDmdtHefsWithoutL(dmdtfs, Tfs, scheme, eqns);
291  addDmdtL(dmdtfs, Tfs, weight, scheme, eqns);
292 }
293 
294 
295 template<class BasePhaseSystem>
297 (
298  const phaseSystem::dmidtfTable& dmidtfs,
299  const phaseSystem::dmdtfTable& Tfs,
300  const latentHeatScheme scheme,
301  phaseSystem::heatTransferTable& eqns
302 ) const
303 {
304  static const dimensionedScalar one(dimless, 1);
305 
306  // Loop the pairs
307  forAllConstIter(phaseSystem::dmidtfTable, dmidtfs, dmidtfIter)
308  {
309  const phasePairKey& key = dmidtfIter.key();
310  const phasePair& pair(this->phasePairs_[key]);
311 
312  const volScalarField& Tf = *Tfs[key];
313 
314  const phaseModel& phase1 = pair.phase1();
315  const phaseModel& phase2 = pair.phase2();
316  const rhoThermo& thermo1 = phase1.thermo();
317  const rhoThermo& thermo2 = phase2.thermo();
318  const basicSpecieMixture* compositionPtr1 =
319  isA<rhoReactionThermo>(thermo1)
320  ? &refCast<const rhoReactionThermo>(thermo1).composition()
321  : static_cast<const basicSpecieMixture*>(nullptr);
322  const basicSpecieMixture* compositionPtr2 =
323  isA<rhoReactionThermo>(thermo2)
324  ? &refCast<const rhoReactionThermo>(thermo2).composition()
325  : static_cast<const basicSpecieMixture*>(nullptr);
326  const volScalarField& he1 = thermo1.he();
327  const volScalarField& he2 = thermo2.he();
328  const volScalarField K1(phase1.K());
329  const volScalarField K2(phase2.K());
330 
331  // Interface enthalpies
332  const volScalarField hsf1(thermo1.hs(thermo1.p(), Tf));
333  const volScalarField hsf2(thermo2.hs(thermo2.p(), Tf));
334 
335  // Loop the species
336  forAllConstIter(HashPtrTable<volScalarField>, *dmidtfIter(), dmidtfJter)
337  {
338  const word& specie = dmidtfJter.key();
339 
340  // Mass transfer rates
341  const volScalarField dmidtf
342  (
343  Pair<word>::compare(pair, key)**dmidtfJter()
344  );
345  const volScalarField dmidtf21(posPart(dmidtf));
346  const volScalarField dmidtf12(negPart(dmidtf));
347 
348  // Specie indices
349  const label speciei1 =
350  compositionPtr1 ? compositionPtr1->species()[specie] : -1;
351  const label speciei2 =
352  compositionPtr2 ? compositionPtr2->species()[specie] : -1;
353 
354  // Interface enthalpies
355  const volScalarField hsfi1
356  (
357  compositionPtr1
358  ? compositionPtr1->Hs(speciei1, thermo1.p(), Tf)
359  : tmp<volScalarField>(hsf1)
360  );
361  const volScalarField hsfi2
362  (
363  compositionPtr2
364  ? compositionPtr2->Hs(speciei2, thermo2.p(), Tf)
365  : tmp<volScalarField>(hsf2)
366  );
367 
368  // Limited mass fractions
369  tmp<volScalarField> tYi1, tYi2;
370  if (this->residualY_ > 0)
371  {
372  tYi1 =
373  compositionPtr1
374  ? max(compositionPtr1->Y(speciei1), this->residualY_)
375  : volScalarField::New("Yi1", this->mesh(), one);
376  tYi2 =
377  compositionPtr2
378  ? max(compositionPtr2->Y(speciei2), this->residualY_)
379  : volScalarField::New("Yi2", this->mesh(), one);
380  }
381 
382  // Transfer of energy from the interface into the bulk
383  switch (scheme)
384  {
385  case latentHeatScheme::symmetric:
386  {
387  *eqns[phase1.name()] += dmidtf*hsfi1;
388  *eqns[phase2.name()] -= dmidtf*hsfi2;
389 
390  break;
391  }
392  case latentHeatScheme::upwind:
393  {
394  // Bulk enthalpies
395  const volScalarField hsi1
396  (
397  compositionPtr1
398  ? compositionPtr1->Hs(speciei1, thermo1.p(), thermo1.T())
399  : thermo1.hs()
400  );
401  const volScalarField hsi2
402  (
403  compositionPtr2
404  ? compositionPtr2->Hs(speciei2, thermo2.p(), thermo2.T())
405  : thermo2.hs()
406  );
407 
408  *eqns[phase1.name()] += dmidtf21*hsfi1 + dmidtf12*hsi1;
409  *eqns[phase2.name()] -= dmidtf12*hsfi2 + dmidtf21*hsi2;
410 
411  break;
412  }
413  }
414  if (this->residualY_ > 0)
415  {
416  *eqns[phase1.name()] +=
417  fvm::Sp(dmidtf12/tYi1(), he1) - dmidtf12/tYi1()*he1;
418  }
419  if (this->residualY_ > 0)
420  {
421  *eqns[phase2.name()] -=
422  fvm::Sp(dmidtf21/tYi2(), he2) - dmidtf21/tYi2()*he2;
423  }
424 
425 
426  // Transfer of kinetic energy
427  *eqns[phase1.name()] += dmidtf21*K2 + dmidtf12*K1;
428  *eqns[phase2.name()] -= dmidtf12*K1 + dmidtf21*K2;
429  }
430  }
431 }
432 
433 
434 template<class BasePhaseSystem>
436 (
437  const phaseSystem::dmidtfTable& dmidtfs,
438  const phaseSystem::dmdtfTable& Tfs,
439  const scalar weight,
440  const latentHeatScheme scheme,
441  phaseSystem::heatTransferTable& eqns
442 ) const
443 {
444  // Loop the pairs
445  forAllConstIter(phaseSystem::dmidtfTable, dmidtfs, dmidtfIter)
446  {
447  const phasePairKey& key = dmidtfIter.key();
448  const phasePair& pair(this->phasePairs_[key]);
449 
450  const volScalarField& Tf = *Tfs[key];
451 
452  const phaseModel& phase1 = pair.phase1();
453  const phaseModel& phase2 = pair.phase2();
454 
455  // Loop the species
456  forAllConstIter(HashPtrTable<volScalarField>, *dmidtfIter(), dmidtfJter)
457  {
458  const word& specie = dmidtfJter.key();
459 
460  // Mass transfer rates
461  const volScalarField dmidtf
462  (
463  Pair<word>::compare(pair, key)**dmidtfJter()
464  );
465  const volScalarField dmidtf21(posPart(dmidtf));
466  const volScalarField dmidtf12(negPart(dmidtf));
467 
468  // Latent heat contribution
469  const volScalarField Li(this->Li(pair, specie, dmidtf, Tf, scheme));
470  *eqns[phase1.name()] +=
471  ((1 - weight)*dmidtf12 + weight*dmidtf21)*Li;
472  *eqns[phase2.name()] +=
473  ((1 - weight)*dmidtf21 + weight*dmidtf12)*Li;
474  }
475  }
476 }
477 
478 
479 template<class BasePhaseSystem>
481 (
482  const phaseSystem::dmidtfTable& dmidtfs,
483  const phaseSystem::dmdtfTable& Tfs,
484  const scalar weight,
485  const latentHeatScheme scheme,
486  phaseSystem::heatTransferTable& eqns
487 ) const
488 {
489  addDmidtHefsWithoutL(dmidtfs, Tfs, scheme, eqns);
490  addDmidtL(dmidtfs, Tfs, weight, scheme, eqns);
491 }
492 
493 
494 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
495 
496 template<class BasePhaseSystem>
498 (
499  const fvMesh& mesh
500 )
501 :
502  heatTransferPhaseSystem(),
503  BasePhaseSystem(mesh),
504  residualY_(this->template lookupOrDefault<scalar>("residualY", -1))
505 {}
506 
507 
508 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
509 
510 template<class BasePhaseSystem>
512 {}
513 
514 
515 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
516 
517 template<class BasePhaseSystem>
520 (
521  const phasePair& pair,
522  const volScalarField& dmdtf,
523  const volScalarField& Tf,
524  const latentHeatScheme scheme
525 ) const
526 {
527  const rhoThermo& thermo1 = pair.phase1().thermo();
528  const rhoThermo& thermo2 = pair.phase2().thermo();
529 
530  // Interface enthalpies
531  const volScalarField haf1(thermo1.ha(thermo1.p(), Tf));
532  const volScalarField haf2(thermo2.ha(thermo2.p(), Tf));
533 
534  switch (scheme)
535  {
536  case latentHeatScheme::symmetric:
537  {
538  return haf2 - haf1;
539  }
540  case latentHeatScheme::upwind:
541  {
542  // Bulk enthalpies
543  const volScalarField ha1(thermo1.ha());
544  const volScalarField ha2(thermo2.ha());
545 
546  return
547  neg0(dmdtf)*haf2 + pos(dmdtf)*ha2
548  - pos0(dmdtf)*haf1 - neg(dmdtf)*ha1;
549  }
550  }
551 
552  return tmp<volScalarField>(nullptr);
553 }
554 
555 
556 template<class BasePhaseSystem>
559 (
560  const phasePair& pair,
561  const scalarField& dmdtf,
562  const scalarField& Tf,
563  const labelUList& cells,
564  const latentHeatScheme scheme
565 ) const
566 {
567  const rhoThermo& thermo1 = pair.phase1().thermo();
568  const rhoThermo& thermo2 = pair.phase2().thermo();
569 
570  // Interface enthalpies
571  const scalarField haf1(thermo1.ha(Tf, cells));
572  const scalarField haf2(thermo2.ha(Tf, cells));
573 
574  switch (scheme)
575  {
576  case latentHeatScheme::symmetric:
577  {
578  return haf2 - haf1;
579  }
580  case latentHeatScheme::upwind:
581  {
582  const scalarField T1(UIndirectList<scalar>(thermo1.T(), cells));
583  const scalarField T2(UIndirectList<scalar>(thermo2.T(), cells));
584 
585  // Bulk enthalpies
586  const scalarField ha1(thermo1.ha(T1, cells));
587  const scalarField ha2(thermo2.ha(T2, cells));
588 
589  return
590  neg0(dmdtf)*haf2 + pos(dmdtf)*ha2
591  - pos0(dmdtf)*haf1 - neg(dmdtf)*ha1;
592  }
593  }
594 
595  return tmp<scalarField>(nullptr);
596 }
597 
598 
599 template<class BasePhaseSystem>
602 (
603  const phasePair& pair,
604  const word& specie,
605  const volScalarField& dmdtf,
606  const volScalarField& Tf,
607  const latentHeatScheme scheme
608 ) const
609 {
610  const rhoThermo& thermo1 = pair.phase1().thermo();
611  const rhoThermo& thermo2 = pair.phase2().thermo();
612  const basicSpecieMixture* compositionPtr1 =
613  isA<rhoReactionThermo>(thermo1)
614  ? &refCast<const rhoReactionThermo>(thermo1).composition()
615  : static_cast<const basicSpecieMixture*>(nullptr);
616  const basicSpecieMixture* compositionPtr2 =
617  isA<rhoReactionThermo>(thermo2)
618  ? &refCast<const rhoReactionThermo>(thermo2).composition()
619  : static_cast<const basicSpecieMixture*>(nullptr);
620  const label speciei1 =
621  compositionPtr1 ? compositionPtr1->species()[specie] : -1;
622  const label speciei2 =
623  compositionPtr2 ? compositionPtr2->species()[specie] : -1;
624 
625  // Interface enthalpies
626  const volScalarField hafi1
627  (
628  compositionPtr1
629  ? compositionPtr1->Ha(speciei1, thermo1.p(), Tf)
630  : thermo1.ha(thermo1.p(), Tf)
631  );
632  const volScalarField hafi2
633  (
634  compositionPtr2
635  ? compositionPtr2->Ha(speciei2, thermo2.p(), Tf)
636  : thermo2.ha(thermo1.p(), Tf)
637  );
638 
639  switch (scheme)
640  {
641  case latentHeatScheme::symmetric:
642  {
643  return hafi2 - hafi1;
644  }
645  case latentHeatScheme::upwind:
646  {
647  // Bulk enthalpies
648  const volScalarField hai1
649  (
650  compositionPtr1
651  ? compositionPtr1->Ha(speciei1, thermo1.p(), thermo1.T())
652  : thermo1.ha()
653  );
654  const volScalarField hai2
655  (
656  compositionPtr2
657  ? compositionPtr2->Ha(speciei2, thermo2.p(), thermo2.T())
658  : thermo2.ha()
659  );
660 
661  return
662  neg0(dmdtf)*hafi2 + pos(dmdtf)*hai2
663  - pos0(dmdtf)*hafi1 - neg(dmdtf)*hai1;
664  }
665  }
666 
667  return tmp<volScalarField>(nullptr);
668 }
669 
670 
671 template<class BasePhaseSystem>
674 (
675  const phasePair& pair,
676  const word& specie,
677  const scalarField& dmdtf,
678  const scalarField& Tf,
679  const labelUList& cells,
680  const latentHeatScheme scheme
681 ) const
682 {
683  const rhoThermo& thermo1 = pair.phase1().thermo();
684  const rhoThermo& thermo2 = pair.phase2().thermo();
685  const basicSpecieMixture* compositionPtr1 =
686  isA<rhoReactionThermo>(thermo1)
687  ? &refCast<const rhoReactionThermo>(thermo1).composition()
688  : static_cast<const basicSpecieMixture*>(nullptr);
689  const basicSpecieMixture* compositionPtr2 =
690  isA<rhoReactionThermo>(thermo2)
691  ? &refCast<const rhoReactionThermo>(thermo2).composition()
692  : static_cast<const basicSpecieMixture*>(nullptr);
693  const label speciei1 =
694  compositionPtr1 ? compositionPtr1->species()[specie] : -1;
695  const label speciei2 =
696  compositionPtr2 ? compositionPtr2->species()[specie] : -1;
697 
698  const scalarField p1(UIndirectList<scalar>(thermo1.p(), cells));
699  const scalarField p2(UIndirectList<scalar>(thermo2.p(), cells));
700 
701  // Interface enthalpies
702  const scalarField hafi1
703  (
704  compositionPtr1
705  ? compositionPtr1->Ha(speciei1, p1, Tf)
706  : thermo1.ha(Tf, cells)
707  );
708  const scalarField hafi2
709  (
710  compositionPtr2
711  ? compositionPtr2->Ha(speciei2, p2, Tf)
712  : thermo2.ha(Tf, cells)
713  );
714 
715  switch (scheme)
716  {
717  case latentHeatScheme::symmetric:
718  {
719  return hafi2 - hafi1;
720  }
721  case latentHeatScheme::upwind:
722  {
723  const scalarField T1(UIndirectList<scalar>(thermo1.T(), cells));
724  const scalarField T2(UIndirectList<scalar>(thermo2.T(), cells));
725 
726  // Bulk enthalpies
727  const scalarField hai1
728  (
729  compositionPtr1
730  ? compositionPtr1->Ha(speciei1, p1, T1)
731  : thermo1.ha(T1, cells)
732  );
733  const scalarField hai2
734  (
735  compositionPtr2
736  ? compositionPtr2->Ha(speciei2, p2, T2)
737  : thermo2.ha(T2, cells)
738  );
739 
740  return
741  neg0(dmdtf)*hafi2 + pos(dmdtf)*hai2
742  - pos0(dmdtf)*hafi1 - neg(dmdtf)*hai1;
743  }
744  }
745 
746  return tmp<scalarField>(nullptr);
747 }
748 
749 
750 template<class BasePhaseSystem>
752 {
753  if (BasePhaseSystem::read())
754  {
755  bool readOK = true;
756 
757  // Models ...
758 
759  return readOK;
760  }
761  else
762  {
763  return false;
764  }
765 }
766 
767 
768 // ************************************************************************* //
virtual tmp< volScalarField > L(const phasePair &pair, const volScalarField &dmdtf, const volScalarField &Tf, const latentHeatScheme scheme) const
Return the latent heat for a given pair, mass transfer rate (used.
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
#define K1
Definition: SHA1.C:167
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define K2
Definition: SHA1.C:168
void addDmdtHefsWithoutL(const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from bulk phase changes,.
static int compare(const Pair< word > &a, const Pair< word > &b)
Compare Pairs.
Definition: Pair.H:153
const dimensionSet dimless
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
mixture thermo1().correctRho(psi1 *(p_rgh - p_rgh_0))
dimensionedScalar posPart(const dimensionedScalar &ds)
dimensionedScalar neg(const dimensionedScalar &ds)
UList< label > labelUList
Definition: UList.H:65
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.
dimensionedScalar pos(const dimensionedScalar &ds)
mixture thermo2().correctRho(psi2 *(p_rgh - p_rgh_0))
dynamicFvMesh & mesh
const cellShapeList & cells
autoPtr< BasicCompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const typename BasicCompressibleMomentumTransportModel::transportModel &transport)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
virtual ~HeatTransferPhaseSystem()
Destructor.
dimensionedScalar neg0(const dimensionedScalar &ds)
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
dimensionedScalar pos0(const dimensionedScalar &ds)
void addDmidtL(const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const scalar weight, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
Add latent heat terms which result from specie phase changes.
void addDmdtL(const phaseSystem::dmdtfTable &dmdtfs, const phaseSystem::dmdtfTable &Tfs, const scalar weight, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
Add latent heat terms which result from bulk phase changes.
const tmp< volScalarField::Internal > & Sp
Definition: alphaSuSp.H:7
void addDmidtHefsWithoutL(const phaseSystem::dmidtfTable &dmidtfs, const phaseSystem::dmdtfTable &Tfs, const latentHeatScheme scheme, phaseSystem::heatTransferTable &eqns) const
Add energy transfer terms which result from specie phase changes,.
virtual tmp< volScalarField > Li(const phasePair &pair, const word &member, const volScalarField &dmidtf, const volScalarField &Tf, const latentHeatScheme scheme) const
Return the latent heat for a given pair, specie, mass transfer rate.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
A class for managing temporary objects.
Definition: PtrList.H:53
HeatTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual bool read()
Read base phaseProperties dictionary.
Calculate the matrix for implicit and explicit sources.
dimensionedScalar negPart(const dimensionedScalar &ds)