40template<
class BasePhaseSystem>
52 if (!
dmdt_.found(iterModel()->pair()))
79template<
class BasePhaseSystem>
98 if (massTransferModels_.found(keyik))
101 massTransferModels_[keyik];
103 word speciesName = interfacePtr->transferSpecie();
105 const word species(speciesName.substr(0, speciesName.find(
'.')));
110 if (massTransferModels_.found(keyki))
113 massTransferModels_[keyki];
115 word speciesName = interfacePtr->transferSpecie();
117 const word species(speciesName.substr(0, speciesName.find(
'.')));
128template<
class BasePhaseSystem>
142 auto& dmdt = tdmdt.ref();
144 if (dmdt_.found(key))
153template<
class BasePhaseSystem>
161 auto& eqn = teqn.ref();
169 for (++iterk; iterk != this->phaseModels_.end(); ++iterk)
171 if (iteri()().
name() != iterk()().
name())
189 auto& dmdtNetki = tdmdtNetki.
ref();
198 auto&
Sp = tSp.ref();
207 auto&
Su = tSu.ref();
210 if (massTransferModels_.found(keyik))
213 massTransferModels_[keyik];
215 dmdtNetki -= *dmdt_[keyik];
241 if (massTransferModels_.found(keyki))
244 massTransferModels_[keyki];
246 dmdtNetki += *dmdt_[keyki];
283template<
class BasePhaseSystem>
291 auto& eqn = teqn.ref();
300 auto&
Sp = tSp.ref();
309 auto&
Su = tSu.ref();
325 if (massTransferModels_.found(key12))
328 massTransferModels_[key12];
338 - this->coeffs(
phase1.name())
339 + this->coeffs(
phase2.name())
351 - this->coeffs(
phase1.name())
352 + this->coeffs(
phase2.name())
362 - this->coeffs(
phase1.name())
363 + this->coeffs(
phase2.name())
375 if (massTransferModels_.found(key21))
378 massTransferModels_[key21];
425template<
class BasePhaseSystem>
437 for (++iterk; iterk != this->phaseModels_.end(); ++iterk)
439 if (iteri()().
name() != iterk()().
name())
449 if (massTransferModels_.found(keyik))
452 massTransferModels_[keyik];
456 *dmdt_[keyik] = Kexp.
ref();
460 if (massTransferModels_.found(keyki))
463 massTransferModels_[keyki];
468 *dmdt_[keyki] = Kexp.
ref();
476template<
class BasePhaseSystem>
485 bool includeDivU(
true);
514 if (massTransferModels_.found(key12))
517 massTransferModels_[key12];
527 includeDivU = interfacePtr->includeDivU();
542 if (massTransferModels_.found(key21))
545 massTransferModels_[key21];
555 includeDivU = interfacePtr->includeDivU();
585 scalar dmdt21 = dmdtNet[celli];
586 scalar coeffs12Cell = coeffs12[celli];
591 SuPhase1[celli] += coeffs1[celli]*dmdt21;
597 if (coeffs12Cell > 0)
600 SpPhase1[celli] -= dmdt21*coeffs12Cell;
602 else if (coeffs12Cell < 0)
606 dmdt21*coeffs12Cell*alpha1Limited;
611 if (coeffs12Cell > 0)
615 dmdt21*coeffs12Cell*alpha1Limited;
617 else if (coeffs12Cell < 0)
620 SpPhase1[celli] -= dmdt21*coeffs12Cell;
628 scalar dmdt12 = -dmdtNet[celli];
629 scalar coeffs21Cell = -coeffs12[celli];
634 SuPhase2[celli] += coeffs2[celli]*dmdt12;
640 if (coeffs21Cell > 0)
643 SpPhase2[celli] -= dmdt12*coeffs21Cell;
645 else if (coeffs21Cell < 0)
649 dmdt12*coeffs21Cell*alpha2Limited;
654 if (coeffs21Cell > 0)
658 coeffs21Cell*dmdt12*alpha2Limited;
660 else if (coeffs21Cell < 0)
663 SpPhase2[celli] -= dmdt12*coeffs21Cell;
672 max(
gMax((dmdt21*coeffs1)()),
gMax((dmdt12*coeffs2)()));
677template<
class BasePhaseSystem>
683 const word speciesName
689 if (iter()->transferSpecie() == speciesName)
700template<
class BasePhaseSystem>
703 bool includeVolChange(
true);
706 if (!iter()->includeVolChange())
708 includeVolChange =
false;
711 return includeVolChange;
const volScalarField & alpha1
const volScalarField & alpha2
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())
DimensionedField< scalar, volMesh > Internal
const GeometricField< Type, PatchField, GeoMesh > & oldTime() const
Return old time field.
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef().
@ NO_REGISTER
Do not request registration (bool: false).
@ REGISTER
Request registration (bool: true).
@ NO_READ
Nothing to be read.
@ AUTO_WRITE
Automatically write from objectRegistry::writeObject().
Defines the attributes of an object for which implicit objectRegistry management is supported,...
static word groupName(StringType base, const word &group)
Create dot-delimited name.group string.
tmp< volScalarField > dmdt(const phasePairKey &key) const
Return total interfacial mass flow rate.
virtual void correctMassSources(const volScalarField &T)
Correct/calculates mass sources dmdt for phases.
MassTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
virtual void massSpeciesTransfer(const Foam::phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName)
Calculate mass transfer for species.
dmdtTable dmdt_
Overall inter-phase mass transfer rates [Kg/s].
HashTable< volScalarField::Internal > SuSpTable
virtual tmp< fvScalarMatrix > volTransfer(const volScalarField &p)
Return the volumetric rate transfer matrix.
virtual bool includeVolChange()
Add volume change in pEq.
massTransferModelTable massTransferModels_
Mass transfer models.
virtual void alphaTransfer(SuSpTable &Su, SuSpTable &Sp)
Calculate mass transfer for alpha's.
virtual tmp< fvScalarMatrix > heatTransfer(const volScalarField &T)
Return the heat transfer matrix.
tmp< volScalarField > calculateL(const volScalarField &dmdtNetki, const phasePairKey &keyik, const phasePairKey &keyki, const volScalarField &T) const
Calculate L between phases.
Pointer management similar to std::unique_ptr, with some additional methods and type checking.
Mesh data needed to do the Finite Volume discretisation.
const word & name() const
The name of this phase.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
const word & name() const
An ordered or unorder pair of phase names. Typically specified as follows.
Description for mass transfer between a pair of phases. The direction of the mass transfer is from th...
const multiphaseInter::phaseModel & phase1() const
const multiphaseInter::phaseModel & phase2() const
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
bool valid() const noexcept
Identical to good(), or bool operator.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
T & ref() const
Return non-const reference to the contents of a non-null managed pointer.
A class for handling words, derived from Foam::string.
Represents 0/1 range or concept. Used for tagged dispatch or clamping.
Fundamental dimensioned constants.
Calculate the divergence of the given field.
Calculate the finiteVolume matrix for implicit and explicit sources.
tmp< GeometricField< Type, fvPatchField, volMesh > > div(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
zeroField Sp(const Foam::zero, const GeometricField< Type, fvPatchField, volMesh > &)
A no-op source.
const dimensionSet dimPressure
dimensionedScalar pos(const dimensionedScalar &ds)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
const dimensionSet dimless
Dimensionless.
const dimensionSet dimEnergy
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
GeometricField< scalar, fvPatchField, volMesh > volScalarField
dimensionSet clamp(const dimensionSet &a, const dimensionSet &range)
GeometricField< scalar, fvsPatchField, surfaceMesh > surfaceScalarField
const dimensionSet dimTemperature(0, 0, 0, 1, 0, 0, 0)
dimensionedScalar neg(const dimensionedScalar &ds)
static constexpr const zero Zero
Global zero (0).
void Sp(fvMatrix< typename Expr::value_type > &m, const Expr2 &mult, const Expr &expression)
const dimensionSet dimVolume(pow3(dimLength))
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimDensity
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
Type gMax(const FieldField< Field, Type > &f)
void Su(fvMatrix< typename Expr::value_type > &m, const Expr &expression)
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
#define forAll(list, i)
Loop across all elements in list.
#define forAllIters(container, iter)
Iterate across all elements in the container object.
#define forAllConstIters(container, iter)
Iterate across all elements of the container object with const access.
const vector L(dict.get< vector >("L"))