Considering the Hertz Knudsen formula, which gives the evaporation-condensation flux based on the kinetic theory for flat interface: More...
Public Member Functions | |
TypeName ("kineticGasEvaporation") | |
kineticGasEvaporation (const dictionary &dict, const phasePair &pair) | |
virtual | ~kineticGasEvaporation ()=default |
virtual tmp< volScalarField > | Kexp (const volScalarField &field) |
virtual tmp< volScalarField > | KSp (label modelVariable, const volScalarField &field) |
virtual tmp< volScalarField > | KSu (label modelVariable, const volScalarField &field) |
virtual const dimensionedScalar & | Tactivate () const noexcept |
virtual bool | includeDivU () const noexcept |
![]() | |
InterfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
~InterfaceCompositionModel ()=default | |
virtual tmp< volScalarField > | dY (const word &speciesName, const volScalarField &Tf) const |
virtual tmp< volScalarField > | Yf (const word &speciesName, const volScalarField &Tf) const |
virtual tmp< volScalarField > | D (const word &speciesName) const |
virtual tmp< volScalarField > | L (const word &speciesName, const volScalarField &Tf) const |
InterfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
~InterfaceCompositionModel ()=default | |
virtual tmp< volScalarField > | dY (const word &speciesName, const volScalarField &Tf) const |
virtual tmp< volScalarField > | D (const word &speciesName) const |
virtual tmp< volScalarField > | L (const word &speciesName, const volScalarField &Tf) const |
virtual void | addMDotL (const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const |
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 |
![]() | |
TypeName ("interfaceCompositionModel") | |
declareRunTimeSelectionTable (autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair)) | |
interfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
virtual | ~interfaceCompositionModel ()=default |
const word | transferSpecie () const |
const phasePair & | pair () const |
bool | includeVolChange () |
const word & | variable () const |
TypeName ("interfaceCompositionModel") | |
declareRunTimeSelectionTable (autoPtr, interfaceCompositionModel, dictionary,(const dictionary &dict, const phasePair &pair),(dict, pair)) | |
interfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
virtual | ~interfaceCompositionModel ()=default |
virtual void | update (const volScalarField &Tf)=0 |
const hashedWordList & | species () const |
bool | transports (word &speciesName) const |
virtual tmp< volScalarField > | YfPrime (const word &speciesName, const volScalarField &Tf) const =0 |
Additional Inherited Members | |
![]() | |
enum | modelVariable { T, P, Y, alpha } |
![]() | |
static autoPtr< interfaceCompositionModel > | New (const dictionary &dict, const phasePair &pair) |
static autoPtr< interfaceCompositionModel > | New (const dictionary &dict, const phasePair &pair) |
![]() | |
template<class ThermoType > | |
const pureMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const pureMixture< ThermoType > &globalThermo) const |
template<class ThermoType > | |
const multiComponentMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const multiComponentMixture< ThermoType > &globalThermo) const |
template<class ThermoType > | |
tmp< volScalarField > | getSpecieMassFraction (const word &speciesName, const pureMixture< ThermoType > &thermo) const |
template<class ThermoType > | |
tmp< volScalarField > | getSpecieMassFraction (const word &speciesName, const multiComponentMixture< ThermoType > &thermo) const |
template<class ThermoType > | |
tmp< volScalarField > | MwMixture (const pureMixture< ThermoType > &thermo) const |
template<class ThermoType > | |
tmp< volScalarField > | MwMixture (const multiComponentMixture< ThermoType > &) const |
template<class ThermoType > | |
const pureMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const pureMixture< ThermoType > &globalThermo) const |
template<class ThermoType > | |
const multiComponentMixture< ThermoType >::thermoType & | getLocalThermo (const word &speciesName, const multiComponentMixture< ThermoType > &globalThermo) const |
![]() | |
const Thermo & | fromThermo_ |
const OtherThermo & | toThermo_ |
const dimensionedScalar | Le_ |
const Thermo & | thermo_ |
const OtherThermo & | otherThermo_ |
![]() | |
modelVariable | modelVariable_ |
bool | includeVolChange_ |
const phasePair & | pair_ |
word | speciesName_ |
const fvMesh & | mesh_ |
const hashedWordList | speciesNames_ |
![]() | |
static const Enum< modelVariable > | modelVariableNames_ |
Considering the Hertz Knudsen formula, which gives the evaporation-condensation flux based on the kinetic theory for flat interface:
where:
![]() | = | mass flux rate [kg/s/m2] |
![]() | = | molecular weight |
![]() | = | saturation temperature |
![]() | = | accommodation coefficient |
![]() | = | universal gas constant |
![]() | = | saturation pressure |
![]() | = | vapor partial pressure |
The Clapeyron-Clausius equation relates the pressure to the temperature for the saturation condition:
where:
![]() | = | latent heat |
![]() | = | inverse of the vapor density |
![]() | = | inverse of the liquid density |
Using the above relations:
This assumes liquid and vapour are in equilibrium, then the accommodation coefficient are equivalent for the interface. This relation is known as the Hertz-Knudsen-Schrage.
Based on the reference:
Example usage:
massTransferModel ( (liquid to gas) { type kineticGasEvaporation; species vapour.gas; C 0.1; isoAlpha 0.1; Tactivate 373; } );
where:
Property | Description | Required | Default value |
---|---|---|---|
C | Coefficient (C > 0 for evaporation, C < 0 for | ||
condensation ) | yes | ||
includeVolChange | Volumen change | no | yes |
isoAlpha | iso-alpha for interface | no | 0.5 |
Tactivate | Saturation temperature | yes | |
species | Specie name on the other phase | no | none |
Definition at line 201 of file kineticGasEvaporation.H.
kineticGasEvaporation | ( | const dictionary & | dict, |
const phasePair & | pair | ||
) |
Definition at line 109 of file kineticGasEvaporation.C.
References Foam::abort(), Foam::constant::electromagnetic::e, Foam::FatalError, FatalErrorInFunction, and IOobject::member().
|
virtualdefault |
TypeName | ( | "kineticGasEvaporation< Thermo, OtherThermo >" | ) |
|
virtual |
Implements interfaceCompositionModel.
Definition at line 186 of file kineticGasEvaporation.C.
References Foam::dimDensity, Foam::dimTemperature, phasePair::from(), L(), Foam::mag(), Foam::max(), IOobject::member(), mesh, Foam::constant::mathematical::pi(), Foam::pow3(), Foam::constant::physicoChemical::R, tmp::ref(), phaseModel::rho(), Foam::sign(), Foam::sqrt(), Foam::T(), T0, phasePair::to(), and Foam::Zero.
|
virtual |
Implements interfaceCompositionModel.
Definition at line 262 of file kineticGasEvaporation.C.
References Foam::pos(), and Foam::sign().
|
virtual |
Implements interfaceCompositionModel.
Definition at line 290 of file kineticGasEvaporation.C.
References Foam::pos(), and Foam::sign().
|
inlinevirtualnoexcept |
Implements interfaceCompositionModel.
Definition at line 278 of file kineticGasEvaporation.H.
|
inlinevirtualnoexcept |
Reimplemented from interfaceCompositionModel.
Definition at line 285 of file kineticGasEvaporation.H.
Copyright © 2011-2018 OpenFOAM | OPENFOAM® is a registered trademark of OpenCFD Ltd.