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(),
67  ),
68  this->mesh(),
70  )
71  );
72  }
73  }
74 }
75 
76 
77 // * * * * * * * * * * * * * * Protected Member Functions * * * * * * * * * //
78 
79 template<class BasePhaseSystem>
82 (
83  const volScalarField& dmdtNetki,
84  const phasePairKey& keyik,
85  const phasePairKey& keyki,
86  const volScalarField& T
87 ) const
88 {
89  auto tL = volScalarField::New
90  (
91  "tL",
93  this->mesh(),
95  );
96  auto& L = tL.ref();
97 
98  if (massTransferModels_.found(keyik))
99  {
100  const autoPtr<interfaceCompositionModel>& interfacePtr =
101  massTransferModels_[keyik];
102 
103  word speciesName = interfacePtr->transferSpecie();
104 
105  const word species(speciesName.substr(0, speciesName.find('.')));
106 
107  L -= neg(dmdtNetki)*interfacePtr->L(species, T);
108  }
109 
110  if (massTransferModels_.found(keyki))
111  {
112  const autoPtr<interfaceCompositionModel>& interfacePtr =
113  massTransferModels_[keyki];
114 
115  word speciesName = interfacePtr->transferSpecie();
116 
117  const word species(speciesName.substr(0, speciesName.find('.')));
118 
119  L -= pos(dmdtNetki)*interfacePtr->L(species, T);
120  }
121 
122  return tL;
123 }
125 
126 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
127 
128 template<class BasePhaseSystem>
131 (
132  const phasePairKey& key
133 ) const
134 {
135  auto tdmdt = volScalarField::New
136  (
137  "dmdt",
139  this->mesh(),
141  );
142  auto& dmdt = tdmdt.ref();
143 
144  if (dmdt_.found(key))
145  {
146  dmdt = *dmdt_[key];
147  }
148 
149  return tdmdt;
150 }
151 
152 
153 template<class BasePhaseSystem>
156 (
157  const volScalarField& T
158 )
159 {
161  auto& eqn = teqn.ref();
162 
163  forAllConstIters(this->phaseModels_, iteri)
164  {
165  const phaseModel& phasei = iteri()();
166 
167  auto iterk = iteri;
168 
169  for (++iterk; iterk != this->phaseModels_.end(); ++iterk)
170  {
171  if (iteri()().name() != iterk()().name())
172  {
173  const phaseModel& phasek = iterk()();
174 
175  // Phase i to phase k
176  const phasePairKey keyik(phasei.name(), phasek.name(), true);
177 
178  // Phase k to phase i
179  const phasePairKey keyki(phasek.name(), phasei.name(), true);
180 
181  // Net mass transfer from k to i phase
182  auto tdmdtNetki = volScalarField::New
183  (
184  "tdmdtYki",
186  this->mesh(),
188  );
189  auto& dmdtNetki = tdmdtNetki.ref();
190 
191  auto tSp = volScalarField::New
192  (
193  "Sp",
195  this->mesh(),
197  );
198  auto& Sp = tSp.ref();
199 
200  auto tSu = volScalarField::New
201  (
202  "Su",
204  this->mesh(),
206  );
207  auto& Su = tSu.ref();
208 
209 
210  if (massTransferModels_.found(keyik))
211  {
212  autoPtr<interfaceCompositionModel>& interfacePtr =
213  massTransferModels_[keyik];
214 
215  dmdtNetki -= *dmdt_[keyik];
216 
217  tmp<volScalarField> KSp =
218  interfacePtr->KSp(interfaceCompositionModel::T, T);
219 
220  if (KSp.valid())
221  {
222  Sp += KSp.ref();
223  }
224 
225  tmp<volScalarField> KSu =
226  interfacePtr->KSu(interfaceCompositionModel::T, T);
227 
228  if (KSu.valid())
229  {
230  Su += KSu.ref();
231  }
232 
233  // If linearization is not provided used full explicit
234  if (!KSp.valid() && !KSu.valid())
235  {
236  Su += *dmdt_[keyik];
237  }
238  }
239 
240  // Looking for mass transfer in the other direction (k to i)
241  if (massTransferModels_.found(keyki))
242  {
243  autoPtr<interfaceCompositionModel>& interfacePtr =
244  massTransferModels_[keyki];
245 
246  dmdtNetki += *dmdt_[keyki];
247 
248 
249  tmp<volScalarField> KSp =
250  interfacePtr->KSp(interfaceCompositionModel::T, T);
251 
252  if (KSp.valid())
253  {
254  Sp += KSp.ref();
255  }
256 
257  tmp<volScalarField> KSu =
258  interfacePtr->KSu(interfaceCompositionModel::T, T);
259 
260  if (KSu.valid())
261  {
262  Su += KSu.ref();
263  }
264 
265  // If linearization is not provided used full explicit
266  if (!KSp.valid() && !KSu.valid())
267  {
268  Su += *dmdt_[keyki];
269  }
270  }
271 
272  tmp<volScalarField> L = calculateL(dmdtNetki, keyik, keyki, T);
273 
274  //eqn -= dmdtNetki*L;
275  eqn -= fvm::Sp(Sp*L.ref(), T) + Su*L.ref();
276  }
277  }
278  }
279  return teqn;
280 }
281 
282 
283 template<class BasePhaseSystem>
286 (
287  const volScalarField& p
288 )
289 {
291  auto& eqn = teqn.ref();
292 
293  auto tSp = volScalarField::New
294  (
295  "Sp",
297  this->mesh(),
299  );
300  auto& Sp = tSp.ref();
301 
302  auto tSu = volScalarField::New
303  (
304  "Su",
306  this->mesh(),
308  );
309  auto& Su = tSu.ref();
310 
311  forAllConstIters(this->totalPhasePairs(), iter)
312  {
313  const phasePair& pair = iter()();
314 
315  const phaseModel& phase1 = pair.phase1();
316  const phaseModel& phase2 = pair.phase2();
317 
318  const phasePairKey key12
319  (
320  phase1.name(),
321  phase2.name(),
322  true
323  );
324 
325  if (massTransferModels_.found(key12))
326  {
327  autoPtr<interfaceCompositionModel>& interfacePtr =
328  massTransferModels_[key12];
329 
330  tmp<volScalarField> KSp =
331  interfacePtr->KSp(interfaceCompositionModel::P, p);
332 
333  if (KSp.valid())
334  {
335  Sp -=
336  KSp.ref()
337  *(
338  - this->coeffs(phase1.name())
339  + this->coeffs(phase2.name())
340  );
341  }
342 
343  tmp<volScalarField> KSu =
344  interfacePtr->KSu(interfaceCompositionModel::P, p);
345 
346  if (KSu.valid())
347  {
348  Su -=
349  KSu.ref()
350  *(
351  - this->coeffs(phase1.name())
352  + this->coeffs(phase2.name())
353  );
354  }
355 
356  // If linearization is not provided used full explicit
357  if (!KSp.valid() && !KSu.valid())
358  {
359  Su -=
360  *dmdt_[key12]
361  *(
362  - this->coeffs(phase1.name())
363  + this->coeffs(phase2.name())
364  );
365  }
366  }
367 
368  const phasePairKey key21
369  (
370  phase2.name(),
371  phase1.name(),
372  true
373  );
374 
375  if (massTransferModels_.found(key21))
376  {
377  autoPtr<interfaceCompositionModel>& interfacePtr =
378  massTransferModels_[key21];
379 
380  tmp<volScalarField> KSp =
381  interfacePtr->KSp(interfaceCompositionModel::P, p);
382 
383  if (KSp.valid())
384  {
385  Sp +=
386  KSp.ref()
387  *(
388  - this->coeffs(phase1.name())
389  + this->coeffs(phase2.name())
390  );
391  }
392 
393  tmp<volScalarField> KSu =
394  interfacePtr->KSu(interfaceCompositionModel::P, p);
395 
396  if (KSu.valid())
397  {
398  Su +=
399  KSu.ref()
400  *(
401  - this->coeffs(phase1.name())
402  + this->coeffs(phase2.name())
403  );
404  }
405 
406  // If linearization is not provided used full explicit
407  if (!KSp.valid() && !KSu.valid())
408  {
409  Su +=
410  *dmdt_[key21]
411  *(
412  - this->coeffs(phase1.name())
413  + this->coeffs(phase2.name())
414  );
415  }
416  }
417 
418  }
419 
420  eqn += fvm::Sp(Sp, p) + Su;
421  return teqn;
422 }
423 
424 
425 template<class BasePhaseSystem>
427 (
428  const volScalarField& T
429 )
430 {
431  forAllConstIters(this->phaseModels_, iteri)
432  {
433  const phaseModel& phasei = iteri()();
434 
435  auto iterk = iteri;
436 
437  for (++iterk; iterk != this->phaseModels_.end(); ++iterk)
438  {
439  if (iteri()().name() != iterk()().name())
440  {
441  const phaseModel& phasek = iterk()();
442 
443  // Phase i to phase k
444  const phasePairKey keyik(phasei.name(), phasek.name(), true);
445 
446  // Phase k to phase i
447  const phasePairKey keyki(phasek.name(), phasei.name(), true);
448 
449  if (massTransferModels_.found(keyik))
450  {
451  autoPtr<interfaceCompositionModel>& interfacePtr =
452  massTransferModels_[keyik];
453 
454  tmp<volScalarField> Kexp = interfacePtr->Kexp(T);
455 
456  *dmdt_[keyik] = Kexp.ref();
457 
458  }
459 
460  if (massTransferModels_.found(keyki))
461  {
462  autoPtr<interfaceCompositionModel>& interfacePtr =
463  massTransferModels_[keyki];
464 
465  // Explicit temperature mass transfer rate
466  const tmp<volScalarField> Kexp = interfacePtr->Kexp(T);
467 
468  *dmdt_[keyki] = Kexp.ref();
469  }
470  }
471  }
472  }
473 }
474 
475 
476 template<class BasePhaseSystem>
478 (
479  SuSpTable& Su,
480  SuSpTable& Sp
481 )
482 {
483  // This term adds/subtracts alpha*div(U) as a source term
484  // for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2)
485  bool includeDivU(true);
486 
487  forAllConstIters(this->totalPhasePairs(), iter)
488  {
489  const phasePair& pair = iter()();
490 
491  const phaseModel& phase1 = pair.phase1();
492  const phaseModel& phase2 = pair.phase2();
493 
494  const volScalarField& alpha1 = pair.phase1();
495  const volScalarField& alpha2 = pair.phase2();
496 
497  tmp<volScalarField> tCoeffs1 = this->coeffs(phase1.name());
498  const volScalarField& coeffs1 = tCoeffs1();
499 
500  tmp<volScalarField> tCoeffs2 = this->coeffs(phase2.name());
501  const volScalarField& coeffs2 = tCoeffs2();
502 
503  // Phase 1 to phase 2
504  const phasePairKey key12
505  (
506  phase1.name(),
507  phase2.name(),
508  true
509  );
510 
511  tmp<volScalarField> tdmdt12(this->dmdt(key12));
512  volScalarField& dmdt12 = tdmdt12.ref();
513 
514  if (massTransferModels_.found(key12))
515  {
516  autoPtr<interfaceCompositionModel>& interfacePtr =
517  massTransferModels_[key12];
518 
519  tmp<volScalarField> KSu =
520  interfacePtr->KSu(interfaceCompositionModel::alpha, phase1);
521 
522  if (KSu.valid())
523  {
524  dmdt12 = KSu.ref();
525  }
526 
527  includeDivU = interfacePtr->includeDivU();
528  }
529 
530 
531  // Phase 2 to phase 1
532  const phasePairKey key21
533  (
534  phase2.name(),
535  phase1.name(),
536  true
537  );
538 
539  tmp<volScalarField> tdmdt21(this->dmdt(key21));
540  volScalarField& dmdt21 = tdmdt21.ref();
541 
542  if (massTransferModels_.found(key21))
543  {
544  autoPtr<interfaceCompositionModel>& interfacePtr =
545  massTransferModels_[key21];
546 
547  tmp<volScalarField> KSu =
548  interfacePtr->KSu(interfaceCompositionModel::alpha, phase2);
549 
550  if (KSu.valid())
551  {
552  dmdt21 = KSu.ref();
553  }
554 
555  includeDivU = interfacePtr->includeDivU();
556  }
557 
558  volScalarField::Internal& SpPhase1 = Sp[phase1.name()];
559 
560  volScalarField::Internal& SuPhase1 = Su[phase1.name()];
561 
562  volScalarField::Internal& SpPhase2 = Sp[phase2.name()];
563 
564  volScalarField::Internal& SuPhase2 = Su[phase2.name()];
565 
566  const volScalarField dmdtNet(dmdt21 - dmdt12);
567 
568  const volScalarField coeffs12(coeffs1 - coeffs2);
569 
570  const surfaceScalarField& phi = this->phi();
571 
572  if (includeDivU)
573  {
574  SuPhase1 += fvc::div(phi)*clamp(alpha1, zero_one{});
575  SuPhase2 += fvc::div(phi)*clamp(alpha2, zero_one{});
576  }
577 
578  // NOTE: dmdtNet is distributed in terms =
579  // Source for phase 1 =
580  // dmdtNet/rho1
581  // - alpha1*dmdtNet(1/rho1 - 1/rho2)
582 
583  forAll(dmdtNet, celli)
584  {
585  scalar dmdt21 = dmdtNet[celli];
586  scalar coeffs12Cell = coeffs12[celli];
587 
588  scalar alpha1Limited = clamp(alpha1[celli], zero_one{});
589 
590  // exp.
591  SuPhase1[celli] += coeffs1[celli]*dmdt21;
592 
593  if (includeDivU)
594  {
595  if (dmdt21 > 0)
596  {
597  if (coeffs12Cell > 0)
598  {
599  // imp
600  SpPhase1[celli] -= dmdt21*coeffs12Cell;
601  }
602  else if (coeffs12Cell < 0)
603  {
604  // exp
605  SuPhase1[celli] -=
606  dmdt21*coeffs12Cell*alpha1Limited;
607  }
608  }
609  else if (dmdt21 < 0)
610  {
611  if (coeffs12Cell > 0)
612  {
613  // exp
614  SuPhase1[celli] -=
615  dmdt21*coeffs12Cell*alpha1Limited;
616  }
617  else if (coeffs12Cell < 0)
618  {
619  // imp
620  SpPhase1[celli] -= dmdt21*coeffs12Cell;
621  }
622  }
623  }
624  }
625 
626  forAll(dmdtNet, celli)
627  {
628  scalar dmdt12 = -dmdtNet[celli];
629  scalar coeffs21Cell = -coeffs12[celli];
630 
631  scalar alpha2Limited = clamp(alpha2[celli], zero_one{});
632 
633  // exp
634  SuPhase2[celli] += coeffs2[celli]*dmdt12;
635 
636  if (includeDivU)
637  {
638  if (dmdt12 > 0)
639  {
640  if (coeffs21Cell > 0)
641  {
642  // imp
643  SpPhase2[celli] -= dmdt12*coeffs21Cell;
644  }
645  else if (coeffs21Cell < 0)
646  {
647  // exp
648  SuPhase2[celli] -=
649  dmdt12*coeffs21Cell*alpha2Limited;
650  }
651  }
652  else if (dmdt12 < 0)
653  {
654  if (coeffs21Cell > 0)
655  {
656  // exp
657  SuPhase2[celli] -=
658  coeffs21Cell*dmdt12*alpha2Limited;
659  }
660  else if (coeffs21Cell < 0)
661  {
662  // imp
663  SpPhase2[celli] -= dmdt12*coeffs21Cell;
664  }
665  }
666  }
667 
668  }
669 
670  // Update ddtAlphaMax
671  this->ddtAlphaMax_ =
672  max(gMax((dmdt21*coeffs1)()), gMax((dmdt12*coeffs2)()));
673  }
674 }
675 
676 
677 template<class BasePhaseSystem>
679 (
680  const phaseModel& phase,
683  const word speciesName
684 )
685 {
686  // Fill the volumetric mass transfer for species
687  forAllConstIters(massTransferModels_, iter)
688  {
689  if (iter()->transferSpecie() == speciesName)
690  {
691  // Explicit source
692  Su +=
693  this->Su()[phase.name()]
694  + this->Sp()[phase.name()]*phase.oldTime();
695  }
696  }
697 }
698 
699 
700 template<class BasePhaseSystem>
702 {
703  bool includeVolChange(true);
704  forAllIters(massTransferModels_, iter)
705  {
706  if (!iter()->includeVolChange())
707  {
708  includeVolChange = false;
709  }
710  }
711  return includeVolChange;
712 }
713 
714 // ************************************************************************* //
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:1370
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
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:72
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.
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
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...
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:180
Request registration (bool: true)
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
Do not request registration (bool: false)
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