Go to the documentation of this file.
37 #ifndef MassTransferPhaseSystem_H
38 #define MassTransferPhaseSystem_H
40 #include "phaseSystem.H"
42 #include "interfaceCompositionModel.H"
53 template<
class BasePhaseSystem>
56 public BasePhaseSystem
150 const word speciesName
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
HashTable< volScalarField::Internal > SuSpTable
Hashing functor for phasePairKey.
A class for handling words, derived from Foam::string.
Single incompressible phase derived from the phase-fraction. Used as part of the multiPhaseMixture fo...
A class for managing temporary objects.
tmp< volScalarField > calculateL(const volScalarField &dmdtNetki, const phasePairKey &keyik, const phasePairKey &keyki, const volScalarField &T) const
Calculate L between phases.
virtual tmp< fvScalarMatrix > volTransfer(const volScalarField &p)
Return the volumetric rate transfer matrix.
auto key(const Type &t) -> typename std::enable_if< std::is_enum< Type >::value, typename std::underlying_type< Type >::type >::type
MassTransferPhaseSystem(const fvMesh &)
Construct from fvMesh.
HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > dmdtTable
virtual void massSpeciesTransfer(const phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName)
Calculate mass transfer for species.
DimensionedField< scalar, volMesh > Internal
Type of the internal field from which this GeometricField is derived.
dmdtTable dmdt_
Overall inter-phase mass transfer rates [Kg/s].
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
virtual ~MassTransferPhaseSystem()=default
Destructor.
GeometricField< scalar, fvPatchField, volMesh > volScalarField
An ordered or unorder pair of phase names. Typically specified as follows.
massTransferModelTable massTransferModels_
Mass transfer models.
Mesh data needed to do the Finite Volume discretisation.
Class for mass transfer between phases.
HashTable< autoPtr< interfaceCompositionModel >, phasePairKey, phasePairKey::hash > massTransferModelTable
virtual tmp< fvScalarMatrix > heatTransfer(const volScalarField &T)
Return the heat transfer matrix.
A HashTable similar to std::unordered_map.
A HashTable of pointers to objects of type <T>, with deallocation management of the pointers.
virtual bool includeVolChange()
Add volume change in pEq.
virtual void correctMassSources(const volScalarField &T)
Correct/calculates mass sources dmdt for phases.
tmp< volScalarField > dmdt(const phasePairKey &key) const
Return total interfacial mass flow rate.
virtual void alphaTransfer(SuSpTable &Su, SuSpTable &Sp)
Calculate mass transfer for alpha's.