Class for mass transfer between phases. More...
#include <MassTransferPhaseSystem.H>


Public Types | |
| typedef HashTable< autoPtr< interfaceCompositionModel >, phasePairKey, phasePairKey::hash > | massTransferModelTable |
| typedef HashTable< volScalarField::Internal > | SuSpTable |
Public Member Functions | |
| MassTransferPhaseSystem (const fvMesh &) | |
| Construct from fvMesh. | |
| virtual | ~MassTransferPhaseSystem ()=default |
| Destructor. | |
| tmp< volScalarField > | dmdt (const phasePairKey &key) const |
| Return total interfacial mass flow rate. | |
| virtual tmp< fvScalarMatrix > | heatTransfer (const volScalarField &T) |
| Return the heat transfer matrix. | |
| virtual tmp< fvScalarMatrix > | volTransfer (const volScalarField &p) |
| Return the volumetric rate transfer matrix. | |
| virtual void | correctMassSources (const volScalarField &T) |
| Correct/calculates mass sources dmdt for phases. | |
| virtual void | alphaTransfer (SuSpTable &Su, SuSpTable &Sp) |
| Calculate mass transfer for alpha's. | |
| virtual void | massSpeciesTransfer (const Foam::phaseModel &phase, volScalarField::Internal &Su, volScalarField::Internal &Sp, const word speciesName) |
| Calculate mass transfer for species. | |
| virtual bool | includeVolChange () |
| Add volume change in pEq. | |
Protected Types | |
| typedef HashPtrTable< volScalarField, phasePairKey, phasePairKey::hash > | dmdtTable |
Protected Member Functions | |
| tmp< volScalarField > | calculateL (const volScalarField &dmdtNetki, const phasePairKey &keyik, const phasePairKey &keyki, const volScalarField &T) const |
| Calculate L between phases. | |
Protected Attributes | |
| dmdtTable | dmdt_ |
| Overall inter-phase mass transfer rates [Kg/s]. | |
| massTransferModelTable | massTransferModels_ |
| Mass transfer models. | |
Class for mass transfer between phases.
Definition at line 51 of file MassTransferPhaseSystem.H.
| typedef HashTable< autoPtr<interfaceCompositionModel>, phasePairKey, phasePairKey::hash > massTransferModelTable |
Definition at line 66 of file MassTransferPhaseSystem.H.
| typedef HashTable<volScalarField::Internal> SuSpTable |
Definition at line 69 of file MassTransferPhaseSystem.H.
|
protected |
Definition at line 82 of file MassTransferPhaseSystem.H.
|
explicit |
Construct from fvMesh.
Definition at line 34 of file MassTransferPhaseSystem.C.
References IOobjectOption::AUTO_WRITE, Foam::dimDensity, Foam::dimTime, dmdt_, forAllConstIters, IOobject::groupName(), massTransferModels_, mesh, name, IOobjectOption::NO_READ, IOobjectOption::REGISTER, timeName, and Foam::Zero.

|
virtualdefault |
Destructor.
|
protected |
Calculate L between phases.
Definition at line 74 of file MassTransferPhaseSystem.C.
References Foam::dimEnergy, Foam::dimMass, L, massTransferModels_, mesh, Foam::neg(), GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, Foam::pos(), Foam::T(), and Foam::Zero.
Referenced by heatTransfer().


| Foam::tmp< Foam::volScalarField > dmdt | ( | const phasePairKey & | key | ) | const |
Return total interfacial mass flow rate.
Definition at line 123 of file MassTransferPhaseSystem.C.
References Foam::dimDensity, Foam::dimTime, dmdt(), dmdt_, mesh, GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, and Foam::Zero.
Referenced by alphaTransfer(), and dmdt().


|
virtual |
Return the heat transfer matrix.
NOTE: Call KSu and KSp with T as variable,if not provided uses dmdt.
Definition at line 148 of file MassTransferPhaseSystem.C.
References calculateL(), Foam::dimDensity, Foam::dimEnergy, Foam::dimTemperature, Foam::dimTime, dmdt_, forAllConstIters, L, massTransferModels_, mesh, phaseModel::name(), name, GeometricField< scalar, fvPatchField, volMesh >::New(), tmp< T >::New(), IOobjectOption::NO_REGISTER, phasei, tmp< T >::ref(), Foam::fvm::Sp(), Foam::Sp(), Foam::Su(), interfaceCompositionModel::T, Foam::T(), tmp< T >::valid(), and Foam::Zero.

|
virtual |
Return the volumetric rate transfer matrix.
NOTE: Call KSu and KSp with p as variable,if not provided uses dmdt.
Definition at line 278 of file MassTransferPhaseSystem.C.
References Foam::dimless, Foam::dimPressure, Foam::dimTime, Foam::dimVolume, dmdt_, forAllConstIters, massTransferModels_, mesh, GeometricField< scalar, fvPatchField, volMesh >::New(), tmp< T >::New(), IOobjectOption::NO_REGISTER, interfaceCompositionModel::P, p, phasePair::phase1(), phase1, phasePair::phase2(), phase2, tmp< T >::ref(), Foam::fvm::Sp(), Foam::Sp(), Foam::Su(), tmp< T >::valid(), and Foam::Zero.

|
virtual |
Correct/calculates mass sources dmdt for phases.
NOTE: Call the kexp() for all the mass transfer models.
Definition at line 419 of file MassTransferPhaseSystem.C.
References dmdt_, forAllConstIters, massTransferModels_, phaseModel::name(), name, phasei, tmp< T >::ref(), and Foam::T().

Calculate mass transfer for alpha's.
Definition at line 470 of file MassTransferPhaseSystem.C.
References interfaceCompositionModel::alpha, alpha1, alpha2, Foam::clamp(), Foam::fvc::div(), dmdt(), forAll, forAllConstIters, Foam::gMax(), massTransferModels_, Foam::max(), phasePair::phase1(), phase1, phasePair::phase2(), phase2, phi, tmp< T >::ref(), Foam::Sp(), Foam::Su(), and tmp< T >::valid().

|
virtual |
Calculate mass transfer for species.
Definition at line 671 of file MassTransferPhaseSystem.C.
References forAllConstIters, massTransferModels_, phaseModel::name(), GeometricField< Type, PatchField, GeoMesh >::oldTime(), Foam::Sp(), Sp, Foam::Su(), and Su.

|
virtual |
Add volume change in pEq.
Definition at line 694 of file MassTransferPhaseSystem.C.
References forAllIters, includeVolChange(), and massTransferModels_.
Referenced by includeVolChange().


|
protected |
Overall inter-phase mass transfer rates [Kg/s].
Definition at line 90 of file MassTransferPhaseSystem.H.
Referenced by correctMassSources(), dmdt(), heatTransfer(), MassTransferPhaseSystem(), and volTransfer().
|
protected |
Mass transfer models.
Definition at line 95 of file MassTransferPhaseSystem.H.
Referenced by alphaTransfer(), calculateL(), correctMassSources(), heatTransfer(), includeVolChange(), massSpeciesTransfer(), MassTransferPhaseSystem(), and volTransfer().