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