Gas diffusion based evaporation/condensation mass transfer model. More...
Public Member Functions | |
TypeName ("diffusionGasEvaporation") | |
Runtime type information. More... | |
diffusionGasEvaporation (const dictionary &dict, const phasePair &pair) | |
Construct from components. More... | |
virtual | ~diffusionGasEvaporation ()=default |
Destructor. More... | |
virtual tmp< volScalarField > | Kexp (const volScalarField &field) |
Explicit total mass transfer coefficient. More... | |
virtual tmp< volScalarField > | KSp (label modelVariable, const volScalarField &field) |
Implicit mass transfer coefficient. More... | |
virtual tmp< volScalarField > | KSu (label modelVariable, const volScalarField &field) |
Explicit mass transfer coefficient. More... | |
virtual const dimensionedScalar & | Tactivate () const noexcept |
Return Tactivate. More... | |
virtual bool | includeDivU () const noexcept |
Add/subtract alpha*div(U) as a source term for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2) More... | |
Public Member Functions inherited from InterfaceCompositionModel< Thermo, OtherThermo > | |
InterfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
Construct from components. More... | |
virtual | ~InterfaceCompositionModel ()=default |
Destructor. More... | |
virtual tmp< volScalarField > | dY (const word &speciesName, const volScalarField &Tf) const |
Mass fraction difference between the interface and the field. More... | |
virtual tmp< volScalarField > | Yf (const word &speciesName, const volScalarField &Tf) const |
Reference mass fraction for species based models. More... | |
virtual tmp< volScalarField > | Dfrom (const word &speciesName) const |
Specie mass diffusivity for pure mixture. More... | |
virtual tmp< volScalarField > | Dto (const word &speciesName) const |
Specie mass diffusivity for specie in a multicomponent. More... | |
virtual tmp< volScalarField > | L (const word &speciesName, const volScalarField &Tf) const |
Latent heat (to - from)(thermo - otherThermo) More... | |
InterfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
Construct from components. More... | |
~InterfaceCompositionModel ()=default | |
Destructor. More... | |
virtual tmp< volScalarField > | dY (const word &speciesName, const volScalarField &Tf) const |
Mass fraction difference between the interface and the field. More... | |
virtual tmp< volScalarField > | D (const word &speciesName) const |
Mass diffusivity. More... | |
virtual tmp< volScalarField > | L (const word &speciesName, const volScalarField &Tf) const |
Latent heat. More... | |
virtual void | addMDotL (const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const |
Add latent heat flow rate to total. More... | |
template<class ThermoType > | |
const Foam::multiComponentMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const multiComponentMixture< ThermoType > &globalThermo) const |
template<class ThermoType > | |
const Foam::pureMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const pureMixture< ThermoType > &globalThermo) const |
template<class ThermoType > | |
Foam::tmp< Foam::volScalarField > | getSpecieMassFraction (const word &speciesName, const multiComponentMixture< ThermoType > &mixture) const |
template<class ThermoType > | |
Foam::tmp< Foam::volScalarField > | getSpecieMassFraction (const word &speciesName, const pureMixture< ThermoType > &mixture) const |
template<class ThermoType > | |
Foam::tmp< Foam::volScalarField > | MwMixture (const pureMixture< ThermoType > &mixture) const |
template<class ThermoType > | |
Foam::tmp< Foam::volScalarField > | MwMixture (const multiComponentMixture< ThermoType > &mixture) const |
Public Member Functions inherited from interfaceCompositionModel | |
TypeName ("interfaceCompositionModel") | |
Runtime type information. More... | |
declareRunTimeSelectionTable (autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair)) | |
interfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
Construct from a dictionary and a phase pair. More... | |
virtual | ~interfaceCompositionModel ()=default |
Destructor. More... | |
const word | transferSpecie () const |
Return the transferring species name. More... | |
const phasePair & | pair () const |
The phase pair. More... | |
const multiphaseInterSystem & | fluid () const |
Return the multiphaseInterSystem this interface belongs to. More... | |
bool | includeVolChange () |
Add volume change in pEq. More... | |
const word & | variable () const |
Returns the variable on which the model is based. More... | |
Additional Inherited Members | |
Public Types inherited from interfaceCompositionModel | |
enum | modelVariable { T, P, Y, alpha } |
Enumeration for variable based mass transfer models. More... | |
Static Public Member Functions inherited from interfaceCompositionModel | |
static autoPtr< interfaceCompositionModel > | New (const dictionary &dict, const phasePair &pair) |
Protected Member Functions inherited from InterfaceCompositionModel< Thermo, OtherThermo > | |
template<class ThermoType > | |
const pureMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const pureMixture< ThermoType > &globalThermo) const |
Get a reference to the local thermo for a pure mixture. More... | |
template<class ThermoType > | |
const multiComponentMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const multiComponentMixture< ThermoType > &globalThermo) const |
Get a reference to the local thermo for a multi component mixture. More... | |
template<class ThermoType > | |
tmp< volScalarField > | getSpecieMassFraction (const word &speciesName, const pureMixture< ThermoType > &thermo) const |
Return mass fraction for a pureMixture equal to one. More... | |
template<class ThermoType > | |
tmp< volScalarField > | getSpecieMassFraction (const word &speciesName, const multiComponentMixture< ThermoType > &thermo) const |
Return mass fraction for speciesName. More... | |
template<class ThermoType > | |
tmp< volScalarField > | MwMixture (const pureMixture< ThermoType > &thermo) const |
Return moleculas weight of the mixture for pureMixture [Kg/mol]. More... | |
template<class ThermoType > | |
tmp< volScalarField > | MwMixture (const multiComponentMixture< ThermoType > &) const |
Return moleculas weight of the mixture for multiComponentMixture. More... | |
template<class ThermoType > | |
const pureMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const pureMixture< ThermoType > &globalThermo) const |
Get a reference to the local thermo for a pure mixture. More... | |
template<class ThermoType > | |
const multiComponentMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const multiComponentMixture< ThermoType > &globalThermo) const |
Get a reference to the local thermo for a multi component mixture. More... | |
Protected Attributes inherited from InterfaceCompositionModel< Thermo, OtherThermo > | |
const Thermo & | fromThermo_ |
Thermo (from) More... | |
const OtherThermo & | toThermo_ |
Other Thermo (to) More... | |
const dimensionedScalar | Le_ |
Lewis number. More... | |
const Thermo & | thermo_ |
Thermo. More... | |
const OtherThermo & | otherThermo_ |
Other Thermo. More... | |
Protected Attributes inherited from interfaceCompositionModel | |
modelVariable | modelVariable_ |
Enumeration for the model variable. More... | |
bool | includeVolChange_ |
Add volume change in pEq. More... | |
const phasePair & | pair_ |
Phase pair. More... | |
word | speciesName_ |
Names of the transferring specie. More... | |
const fvMesh & | mesh_ |
Reference to mesh. More... | |
Static Protected Attributes inherited from interfaceCompositionModel | |
static const Enum< modelVariable > | modelVariableNames_ |
Selection names for the modelVariable. More... | |
Gas diffusion based evaporation/condensation mass transfer model.
THE evaporation rate is given by:
where:
= | mass flux rate [kg/s/m2] | |
= | gas phase density | |
= | diffusion coefficient | |
= | model coefficient | |
= | normal derivative of evaporated component | |
= | mass fraction at the surface |
massTransferModel ( (liquid to gas) { // Mandatory entries type diffusionGasEvaporation; species vapour.gas; C 1; saturationPressure { type Antoine; A 3.55; B 643; C -198; } // Optional entries isoAlpha <scalar>; Tactivate <scalar>; // Inherited entries ... } );
where:
Property | Description | Type | Reqd | Deflt |
---|---|---|---|---|
type | Type name: diffusionGasEvaporation | word | yes | - |
saturationPressure | Saturation model | dict | yes | - |
isoAlpha | Iso-alpha for interface | scalar | no | 0.5 |
C | Model coefficient | scalar | yes | - |
Tactivate | Saturation temperature | scalar | no | 0 |
The inherited entries are elaborated in:
Definition at line 181 of file diffusionGasEvaporation.H.
diffusionGasEvaporation | ( | const dictionary & | dict, |
const phasePair & | pair | ||
) |
Construct from components.
Definition at line 66 of file diffusionGasEvaporation.C.
|
virtualdefault |
Destructor.
TypeName | ( | "diffusionGasEvaporation< Thermo, OtherThermo >" | ) |
Runtime type information.
|
virtual |
Explicit total mass transfer coefficient.
Implements interfaceCompositionModel.
Definition at line 118 of file diffusionGasEvaporation.C.
References Foam::ensightOutput::debug, Foam::dimDensity, Foam::dimLength, Foam::dimMass, Foam::dimMoles, Foam::dimTime, fluid, phasePair::from(), Foam::fvc::grad(), Foam::max(), IOobject::member(), mesh, tmp< T >::New(), GeometricField< scalar, fvPatchField, volMesh >::New(), IOobjectOption::NO_REGISTER, Foam::pos(), phaseModel::rho(), rhog(), Foam::sqr(), Foam::T(), phasePair::to(), regIOobject::write(), and Foam::Zero.
|
virtual |
Implicit mass transfer coefficient.
Implements interfaceCompositionModel.
Definition at line 214 of file diffusionGasEvaporation.C.
|
virtual |
Explicit mass transfer coefficient.
Implements interfaceCompositionModel.
Definition at line 227 of file diffusionGasEvaporation.C.
|
inlinevirtualnoexcept |
Return Tactivate.
Implements interfaceCompositionModel.
Definition at line 283 of file diffusionGasEvaporation.H.
|
inlinevirtualnoexcept |
Add/subtract alpha*div(U) as a source term for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2)
Reimplemented from interfaceCompositionModel.
Definition at line 292 of file diffusionGasEvaporation.H.