MassTransferPhaseSystem.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2017-2022 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
11  This file is part of OpenFOAM.
12 
13  OpenFOAM is free software: you can redistribute it and/or modify it
14  under the terms of the GNU General Public License as published by
15  the Free Software Foundation, either version 3 of the License, or
16  (at your option) any later version.
17 
18  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21  for more details.
22 
23  You should have received a copy of the GNU General Public License
24  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25 
26 \*---------------------------------------------------------------------------*/
27 
29 #include "HashPtrTable.H"
30 #include "fvcDiv.H"
31 #include "fvmSup.H"
32 #include "fvMatrix.H"
33 #include "volFields.H"
34 #include "fundamentalConstants.H"
35 
36 using namespace Foam::multiphaseInter;
37 
38 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
39 
40 template<class BasePhaseSystem>
42 (
43  const fvMesh& mesh
44 )
45 :
46  BasePhaseSystem(mesh)
47 {
48  this->generatePairsAndSubModels("massTransferModel", massTransferModels_);
49 
51  {
52  if (!dmdt_.found(iterModel()->pair()))
53  {
54  dmdt_.set
55  (
56  iterModel()->pair(),
57  new volScalarField
58  (
59  IOobject
60  (
61  IOobject::groupName("dmdt",iterModel()->pair().name()),
62  this->mesh().time().timeName(),
63  this->mesh(),
66  ),
67  this->mesh(),
69  )
70  );
71  }
72  }
73 }
74 
75 
76 // * * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * //
77 
78 template<class BasePhaseSystem>
81 (
82  const volScalarField& dmdtNetki,
83  const phasePairKey& keyik,
84  const phasePairKey& keyki,
85  const volScalarField& T
86 ) const
87 {
88  auto tL = tmp<volScalarField>::New
89  (
90  IOobject
91  (
92  "tL",
93  this->mesh().time().timeName(),
94  this->mesh(),
97  ),
98  this->mesh(),
100  );
101  auto& L = tL.ref();
102 
103  if (massTransferModels_.found(keyik))
104  {
105  const autoPtr<interfaceCompositionModel>& interfacePtr =
106  massTransferModels_[keyik];
107 
108  word speciesName = interfacePtr->transferSpecie();
109 
110  const word species(speciesName.substr(0, speciesName.find('.')));
111 
112  L -= neg(dmdtNetki)*interfacePtr->L(species, T);
113  }
114 
115  if (massTransferModels_.found(keyki))
116  {
117  const autoPtr<interfaceCompositionModel>& interfacePtr =
118  massTransferModels_[keyki];
119 
120  word speciesName = interfacePtr->transferSpecie();
121 
122  const word species(speciesName.substr(0, speciesName.find('.')));
123 
124  L -= pos(dmdtNetki)*interfacePtr->L(species, T);
125  }
126 
127  return tL;
128 }
130 
131 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
132 
133 template<class BasePhaseSystem>
136 (
137  const phasePairKey& key
138 ) const
139 {
140  auto tdmdt = tmp<volScalarField>::New
141  (
142  IOobject
143  (
144  "dmdt",
145  this->mesh().time().timeName(),
146  this->mesh()
147  ),
148  this->mesh(),
150  );
151  auto& dmdt = tdmdt.ref();
152 
153  if (dmdt_.found(key))
154  {
155  dmdt = *dmdt_[key];
156  }
157 
158  return tdmdt;
159 }
160 
161 
162 template<class BasePhaseSystem>
165 (
166  const volScalarField& T
167 )
168 {
170  auto& eqn = teqn.ref();
171 
172  forAllConstIters(this->phaseModels_, iteri)
173  {
174  const phaseModel& phasei = iteri()();
175 
176  auto iterk = iteri;
177 
178  for (++iterk; iterk != this->phaseModels_.end(); ++iterk)
179  {
180  if (iteri()().name() != iterk()().name())
181  {
182  const phaseModel& phasek = iterk()();
183 
184  // Phase i to phase k
185  const phasePairKey keyik(phasei.name(), phasek.name(), true);
186 
187  // Phase k to phase i
188  const phasePairKey keyki(phasek.name(), phasei.name(), true);
189 
190  // Net mass transfer from k to i phase
191  auto tdmdtNetki = tmp<volScalarField>::New
192  (
193  IOobject
194  (
195  "tdmdtYki",
196  this->mesh().time().timeName(),
197  this->mesh()
198  ),
199  this->mesh(),
201  );
202  auto& dmdtNetki = tdmdtNetki.ref();
203 
204  auto tSp = tmp<volScalarField>::New
205  (
206  IOobject
207  (
208  "Sp",
209  this->mesh().time().timeName(),
210  this->mesh()
211  ),
212  this->mesh(),
214  );
215  auto& Sp = tSp.ref();
216 
217  auto tSu = tmp<volScalarField>::New
218  (
219  IOobject
220  (
221  "Su",
222  this->mesh().time().timeName(),
223  this->mesh()
224  ),
225  this->mesh(),
227  );
228  auto& Su = tSu.ref();
229 
230 
231  if (massTransferModels_.found(keyik))
232  {
233  autoPtr<interfaceCompositionModel>& interfacePtr =
234  massTransferModels_[keyik];
235 
236  dmdtNetki -= *dmdt_[keyik];
237 
238  tmp<volScalarField> KSp =
239  interfacePtr->KSp(interfaceCompositionModel::T, T);
240 
241  if (KSp.valid())
242  {
243  Sp += KSp.ref();
244  }
245 
246  tmp<volScalarField> KSu =
247  interfacePtr->KSu(interfaceCompositionModel::T, T);
248 
249  if (KSu.valid())
250  {
251  Su += KSu.ref();
252  }
253 
254  // If linearization is not provided used full explicit
255  if (!KSp.valid() && !KSu.valid())
256  {
257  Su += *dmdt_[keyik];
258  }
259  }
260 
261  // Looking for mass transfer in the other direction (k to i)
262  if (massTransferModels_.found(keyki))
263  {
264  autoPtr<interfaceCompositionModel>& interfacePtr =
265  massTransferModels_[keyki];
266 
267  dmdtNetki += *dmdt_[keyki];
268 
269 
270  tmp<volScalarField> KSp =
271  interfacePtr->KSp(interfaceCompositionModel::T, T);
272 
273  if (KSp.valid())
274  {
275  Sp += KSp.ref();
276  }
277 
278  tmp<volScalarField> KSu =
279  interfacePtr->KSu(interfaceCompositionModel::T, T);
280 
281  if (KSu.valid())
282  {
283  Su += KSu.ref();
284  }
285 
286  // If linearization is not provided used full explicit
287  if (!KSp.valid() && !KSu.valid())
288  {
289  Su += *dmdt_[keyki];
290  }
291  }
292 
293  tmp<volScalarField> L = calculateL(dmdtNetki, keyik, keyki, T);
294 
295  //eqn -= dmdtNetki*L;
296  eqn -= fvm::Sp(Sp*L.ref(), T) + Su*L.ref();
297  }
298  }
299  }
300  return teqn;
301 }
302 
303 
304 template<class BasePhaseSystem>
307 (
308  const volScalarField& p
309 )
310 {
312  auto& eqn = teqn.ref();
313 
314  auto tSp = tmp<volScalarField>::New
315  (
316  IOobject
317  (
318  "Sp",
319  this->mesh().time().timeName(),
320  this->mesh()
321  ),
322  this->mesh(),
324  );
325  auto& Sp = tSp.ref();
326 
327  auto tSu = tmp<volScalarField>::New
328  (
329  IOobject
330  (
331  "Su",
332  this->mesh().time().timeName(),
333  this->mesh()
334  ),
335  this->mesh(),
337  );
338  auto& Su = tSu.ref();
339 
340  forAllConstIters(this->totalPhasePairs(), iter)
341  {
342  const phasePair& pair = iter()();
343 
344  const phaseModel& phase1 = pair.phase1();
345  const phaseModel& phase2 = pair.phase2();
346 
347  const phasePairKey key12
348  (
349  phase1.name(),
350  phase2.name(),
351  true
352  );
353 
354  if (massTransferModels_.found(key12))
355  {
356  autoPtr<interfaceCompositionModel>& interfacePtr =
357  massTransferModels_[key12];
358 
359  tmp<volScalarField> KSp =
360  interfacePtr->KSp(interfaceCompositionModel::P, p);
361 
362  if (KSp.valid())
363  {
364  Sp -=
365  KSp.ref()
366  *(
367  - this->coeffs(phase1.name())
368  + this->coeffs(phase2.name())
369  );
370  }
371 
372  tmp<volScalarField> KSu =
373  interfacePtr->KSu(interfaceCompositionModel::P, p);
374 
375  if (KSu.valid())
376  {
377  Su -=
378  KSu.ref()
379  *(
380  - this->coeffs(phase1.name())
381  + this->coeffs(phase2.name())
382  );
383  }
384 
385  // If linearization is not provided used full explicit
386  if (!KSp.valid() && !KSu.valid())
387  {
388  Su -=
389  *dmdt_[key12]
390  *(
391  - this->coeffs(phase1.name())
392  + this->coeffs(phase2.name())
393  );
394  }
395  }
396 
397  const phasePairKey key21
398  (
399  phase2.name(),
400  phase1.name(),
401  true
402  );
403 
404  if (massTransferModels_.found(key21))
405  {
406  autoPtr<interfaceCompositionModel>& interfacePtr =
407  massTransferModels_[key21];
408 
409  tmp<volScalarField> KSp =
410  interfacePtr->KSp(interfaceCompositionModel::P, p);
411 
412  if (KSp.valid())
413  {
414  Sp +=
415  KSp.ref()
416  *(
417  - this->coeffs(phase1.name())
418  + this->coeffs(phase2.name())
419  );
420  }
421 
422  tmp<volScalarField> KSu =
423  interfacePtr->KSu(interfaceCompositionModel::P, p);
424 
425  if (KSu.valid())
426  {
427  Su +=
428  KSu.ref()
429  *(
430  - this->coeffs(phase1.name())
431  + this->coeffs(phase2.name())
432  );
433  }
434 
435  // If linearization is not provided used full explicit
436  if (!KSp.valid() && !KSu.valid())
437  {
438  Su +=
439  *dmdt_[key21]
440  *(
441  - this->coeffs(phase1.name())
442  + this->coeffs(phase2.name())
443  );
444  }
445  }
446 
447  }
448 
449  eqn += fvm::Sp(Sp, p) + Su;
450  return teqn;
451 }
452 
453 
454 template<class BasePhaseSystem>
456 (
457  const volScalarField& T
458 )
459 {
460  forAllConstIters(this->phaseModels_, iteri)
461  {
462  const phaseModel& phasei = iteri()();
463 
464  auto iterk = iteri;
465 
466  for (++iterk; iterk != this->phaseModels_.end(); ++iterk)
467  {
468  if (iteri()().name() != iterk()().name())
469  {
470  const phaseModel& phasek = iterk()();
471 
472  // Phase i to phase k
473  const phasePairKey keyik(phasei.name(), phasek.name(), true);
474 
475  // Phase k to phase i
476  const phasePairKey keyki(phasek.name(), phasei.name(), true);
477 
478  if (massTransferModels_.found(keyik))
479  {
480  autoPtr<interfaceCompositionModel>& interfacePtr =
481  massTransferModels_[keyik];
482 
483  tmp<volScalarField> Kexp = interfacePtr->Kexp(T);
484 
485  *dmdt_[keyik] = Kexp.ref();
486 
487  }
488 
489  if (massTransferModels_.found(keyki))
490  {
491  autoPtr<interfaceCompositionModel>& interfacePtr =
492  massTransferModels_[keyki];
493 
494  // Explicit temperature mass transfer rate
495  const tmp<volScalarField> Kexp = interfacePtr->Kexp(T);
496 
497  *dmdt_[keyki] = Kexp.ref();
498  }
499  }
500  }
501  }
502 }
503 
504 
505 template<class BasePhaseSystem>
507 (
508  SuSpTable& Su,
509  SuSpTable& Sp
510 )
511 {
512  // This term adds/subtracts alpha*div(U) as a source term
513  // for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2)
514  bool includeDivU(true);
515 
516  forAllConstIters(this->totalPhasePairs(), iter)
517  {
518  const phasePair& pair = iter()();
519 
520  const phaseModel& phase1 = pair.phase1();
521  const phaseModel& phase2 = pair.phase2();
522 
523  const volScalarField& alpha1 = pair.phase1();
524  const volScalarField& alpha2 = pair.phase2();
525 
526  tmp<volScalarField> tCoeffs1 = this->coeffs(phase1.name());
527  const volScalarField& coeffs1 = tCoeffs1();
528 
529  tmp<volScalarField> tCoeffs2 = this->coeffs(phase2.name());
530  const volScalarField& coeffs2 = tCoeffs2();
531 
532  // Phase 1 to phase 2
533  const phasePairKey key12
534  (
535  phase1.name(),
536  phase2.name(),
537  true
538  );
539 
540  tmp<volScalarField> tdmdt12(this->dmdt(key12));
541  volScalarField& dmdt12 = tdmdt12.ref();
542 
543  if (massTransferModels_.found(key12))
544  {
545  autoPtr<interfaceCompositionModel>& interfacePtr =
546  massTransferModels_[key12];
547 
548  tmp<volScalarField> KSu =
549  interfacePtr->KSu(interfaceCompositionModel::alpha, phase1);
550 
551  if (KSu.valid())
552  {
553  dmdt12 = KSu.ref();
554  }
555 
556  includeDivU = interfacePtr->includeDivU();
557  }
558 
559 
560  // Phase 2 to phase 1
561  const phasePairKey key21
562  (
563  phase2.name(),
564  phase1.name(),
565  true
566  );
567 
568  tmp<volScalarField> tdmdt21(this->dmdt(key21));
569  volScalarField& dmdt21 = tdmdt21.ref();
570 
571  if (massTransferModels_.found(key21))
572  {
573  autoPtr<interfaceCompositionModel>& interfacePtr =
574  massTransferModels_[key21];
575 
576  tmp<volScalarField> KSu =
577  interfacePtr->KSu(interfaceCompositionModel::alpha, phase2);
578 
579  if (KSu.valid())
580  {
581  dmdt21 = KSu.ref();
582  }
583 
584  includeDivU = interfacePtr->includeDivU();
585  }
586 
587  volScalarField::Internal& SpPhase1 = Sp[phase1.name()];
588 
589  volScalarField::Internal& SuPhase1 = Su[phase1.name()];
590 
591  volScalarField::Internal& SpPhase2 = Sp[phase2.name()];
592 
593  volScalarField::Internal& SuPhase2 = Su[phase2.name()];
594 
595  const volScalarField dmdtNet(dmdt21 - dmdt12);
596 
597  const volScalarField coeffs12(coeffs1 - coeffs2);
598 
599  const surfaceScalarField& phi = this->phi();
600 
601  if (includeDivU)
602  {
603  SuPhase1 += fvc::div(phi)*clamp(alpha1, zero_one{});
604  SuPhase2 += fvc::div(phi)*clamp(alpha2, zero_one{});
605  }
606 
607  // NOTE: dmdtNet is distributed in terms =
608  // Source for phase 1 =
609  // dmdtNet/rho1
610  // - alpha1*dmdtNet(1/rho1 - 1/rho2)
611 
612  forAll(dmdtNet, celli)
613  {
614  scalar dmdt21 = dmdtNet[celli];
615  scalar coeffs12Cell = coeffs12[celli];
616 
617  scalar alpha1Limited = clamp(alpha1[celli], zero_one{});
618 
619  // exp.
620  SuPhase1[celli] += coeffs1[celli]*dmdt21;
621 
622  if (includeDivU)
623  {
624  if (dmdt21 > 0)
625  {
626  if (coeffs12Cell > 0)
627  {
628  // imp
629  SpPhase1[celli] -= dmdt21*coeffs12Cell;
630  }
631  else if (coeffs12Cell < 0)
632  {
633  // exp
634  SuPhase1[celli] -=
635  dmdt21*coeffs12Cell*alpha1Limited;
636  }
637  }
638  else if (dmdt21 < 0)
639  {
640  if (coeffs12Cell > 0)
641  {
642  // exp
643  SuPhase1[celli] -=
644  dmdt21*coeffs12Cell*alpha1Limited;
645  }
646  else if (coeffs12Cell < 0)
647  {
648  // imp
649  SpPhase1[celli] -= dmdt21*coeffs12Cell;
650  }
651  }
652  }
653  }
654 
655  forAll(dmdtNet, celli)
656  {
657  scalar dmdt12 = -dmdtNet[celli];
658  scalar coeffs21Cell = -coeffs12[celli];
659 
660  scalar alpha2Limited = clamp(alpha2[celli], zero_one{});
661 
662  // exp
663  SuPhase2[celli] += coeffs2[celli]*dmdt12;
664 
665  if (includeDivU)
666  {
667  if (dmdt12 > 0)
668  {
669  if (coeffs21Cell > 0)
670  {
671  // imp
672  SpPhase2[celli] -= dmdt12*coeffs21Cell;
673  }
674  else if (coeffs21Cell < 0)
675  {
676  // exp
677  SuPhase2[celli] -=
678  dmdt12*coeffs21Cell*alpha2Limited;
679  }
680  }
681  else if (dmdt12 < 0)
682  {
683  if (coeffs21Cell > 0)
684  {
685  // exp
686  SuPhase2[celli] -=
687  coeffs21Cell*dmdt12*alpha2Limited;
688  }
689  else if (coeffs21Cell < 0)
690  {
691  // imp
692  SpPhase2[celli] -= dmdt12*coeffs21Cell;
693  }
694  }
695  }
696 
697  }
698 
699  // Update ddtAlphaMax
700  this->ddtAlphaMax_ =
701  max(gMax((dmdt21*coeffs1)()), gMax((dmdt12*coeffs2)()));
702  }
703 }
704 
705 
706 template<class BasePhaseSystem>
708 (
709  const phaseModel& phase,
712  const word speciesName
713 )
714 {
715  // Fill the volumetric mass transfer for species
716  forAllConstIters(massTransferModels_, iter)
717  {
718  if (iter()->transferSpecie() == speciesName)
719  {
720  // Explicit source
721  Su +=
722  this->Su()[phase.name()]
723  + this->Sp()[phase.name()]*phase.oldTime();
724  }
725  }
726 }
727 
728 
729 template<class BasePhaseSystem>
731 {
732  bool includeVolChange(true);
733  forAllIters(massTransferModels_, iter)
734  {
735  if (!iter()->includeVolChange())
736  {
737  includeVolChange = false;
738  }
739  }
740  return includeVolChange;
741 }
742 
743 // ************************************************************************* //
bool valid() const noexcept
Identical to good(), or bool operator.
Definition: tmp.H:472
virtual bool includeVolChange()
Add volume change in pEq.
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
const multiphaseInter::phaseModel & phase2() const
Definition: phasePairI.H:29
zeroField Su
Definition: alphaSuSp.H:1
label phasei
Definition: pEqn.H:27
const vector L(dict.get< vector >("L"))
virtual tmp< volScalarField > KSp(label modelVariable, const volScalarField &field)=0
Implicit mass transfer.
bool found(const Key &key) const
Same as contains()
Definition: HashTable.H:1354
const multiphaseInter::phaseModel & phase1() const
Definition: phasePairI.H:23
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
Definition: tmpI.H:235
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
virtual tmp< fvScalarMatrix > heatTransfer(const volScalarField &T)
Return the heat transfer matrix.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcDiv.C:42
dmdtTable dmdt_
Overall inter-phase mass transfer rates [Kg/s].
virtual tmp< volScalarField > KSu(label modelVariable, const volScalarField &field)=0
Explicit mass transfer.
const word & name() const
The name of this phase.
Definition: phaseModel.H:122
const volScalarField & alpha2
Ignore writing from objectRegistry::writeObject()
const dimensionSet dimless
Dimensionless.
phaseModel & phase2
tmp< volScalarField > dmdt(const phasePairKey &key) const
Return total interfacial mass flow rate.
dimensionedScalar neg(const dimensionedScalar &ds)
Fundamental dimensioned constants.
virtual void correctMassSources(const volScalarField &T)
Correct/calculates mass sources dmdt for phases.
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
GeometricField< scalar, fvPatchField, volMesh > volScalarField
Definition: volFieldsFwd.H:81
word timeName
Definition: getTimeIndex.H:3
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
const word transferSpecie() const
Return the transferring species name.
virtual tmp< fvScalarMatrix > volTransfer(const volScalarField &p)
Return the volumetric rate transfer matrix.
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
A class for handling words, derived from Foam::string.
Definition: word.H:63
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
Definition: dimensionSets.H:52
bool set(const Key &key, T *ptr)
Assign a new entry, overwrites existing.
virtual tmp< volScalarField > Kexp(const volScalarField &field)=0
Explicit full mass transfer.
#define forAllIters(container, iter)
Iterate across all elements in the container object.
Definition: stdFoam.H:336
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
const dimensionSet dimPressure
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
Calculate the divergence of the given field.
An ordered or unorder pair of phase names. Typically specified as follows.
Definition: phasePairKey.H:61
phaseModel & phase1
Type gMax(const FieldField< Field, Type > &f)
MassTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
const dimensionSet dimEnergy
const dimensionSet dimDensity
virtual void alphaTransfer(SuSpTable &Su, SuSpTable &Sp)
Calculate mass transfer for alpha&#39;s.
virtual void massSpeciesTransfer(const Foam::phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName)
Calculate mass transfer for species.
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
Definition: foamGltfBase.H:103
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: areaFieldsFwd.H:42
Nothing to be read.
Automatically write from objectRegistry::writeObject()
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phase.H:50
const word & name() const
Definition: phaseModel.H:166
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
Pointer management similar to std::unique_ptr, with some additional methods and type checking...
Definition: HashPtrTable.H:48
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
volScalarField & p
A class for managing temporary objects.
Definition: HashPtrTable.H:50
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
virtual tmp< volScalarField > L(const word &speciesName, const volScalarField &Tf) const =0
Latent heat (delta Hc)
tmp< volScalarField > calculateL(const volScalarField &dmdtNetki, const phasePairKey &keyik, const phasePairKey &keyki, const volScalarField &T) const
Calculate L between phases.
Defines the attributes of an object for which implicit objectRegistry management is supported...
Definition: IOobject.H:172
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
Definition: phaseModel.H:53
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
Definition: phasePair.H:49
Calculate the finiteVolume matrix for implicit and explicit sources.
forAllConstIters(mixture.phases(), phase)
Definition: pEqn.H:28
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
Definition: dimensionSet.C:271
zeroField Sp
Definition: alphaSuSp.H:2
massTransferModelTable massTransferModels_
Mass transfer models.
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127
const volScalarField & alpha1