kineticGasEvaporation< Thermo, OtherThermo > Class Template Reference

Considering the Hertz Knudsen formula, which gives the evaporation-condensation flux based on the kinetic theory for flat interface: More...

Inheritance diagram for kineticGasEvaporation< Thermo, OtherThermo >:
[legend]
Collaboration diagram for kineticGasEvaporation< Thermo, OtherThermo >:
[legend]

Public Member Functions

 TypeName ("kineticGasEvaporation")
 Runtime type information. More...
 
 kineticGasEvaporation (const dictionary &dict, const phasePair &pair)
 Construct from components. More...
 
virtual ~kineticGasEvaporation ()=default
 Destructor. More...
 
virtual tmp< volScalarFieldKexp (const volScalarField &field)
 Explicit total mass transfer coefficient. More...
 
virtual tmp< volScalarFieldKSp (label modelVariable, const volScalarField &field)
 Implicit mass transfer coefficient. More...
 
virtual tmp< volScalarFieldKSu (label modelVariable, const volScalarField &field)
 Explicit mass transfer coefficient. More...
 
virtual const dimensionedScalarTactivate () const noexcept
 Return Tactivate. More...
 
virtual bool includeDivU () const noexcept
 
- 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< volScalarFielddY (const word &speciesName, const volScalarField &Tf) const
 Mass fraction difference between the interface and the field. More...
 
virtual tmp< volScalarFieldYf (const word &speciesName, const volScalarField &Tf) const
 Reference mass fraction for species based models. More...
 
virtual tmp< volScalarFieldDfrom (const word &speciesName) const
 Specie mass diffusivity for pure mixture. More...
 
virtual tmp< volScalarFieldDto (const word &speciesName) const
 Specie mass diffusivity for specie in a multicomponent. More...
 
virtual tmp< volScalarFieldL (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< volScalarFielddY (const word &speciesName, const volScalarField &Tf) const
 Mass fraction difference between the interface and the field. More...
 
virtual tmp< volScalarFieldD (const word &speciesName) const
 Mass diffusivity. More...
 
virtual tmp< volScalarFieldL (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::volScalarFieldgetSpecieMassFraction (const word &speciesName, const multiComponentMixture< ThermoType > &mixture) const
 
template<class ThermoType >
Foam::tmp< Foam::volScalarFieldgetSpecieMassFraction (const word &speciesName, const pureMixture< ThermoType > &mixture) const
 
template<class ThermoType >
Foam::tmp< Foam::volScalarFieldMwMixture (const pureMixture< ThermoType > &mixture) const
 
template<class ThermoType >
Foam::tmp< Foam::volScalarFieldMwMixture (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 phasePairpair () const
 The phase pair. More...
 
const multiphaseInterSystemfluid () const
 Return the multiphaseInterSystem this interface belongs to. More...
 
virtual tmp< volScalarFieldYf (const word &speciesName, const volScalarField &Tf) const =0
 Interface mass fraction. More...
 
virtual tmp< volScalarFielddY (const word &speciesName, const volScalarField &Tf) const =0
 Mass fraction difference between the interface and the field. More...
 
virtual tmp< volScalarFieldDfrom (const word &speciesName) const =0
 Specie mass diffusivity for pure mixture. More...
 
virtual tmp< volScalarFieldDto (const word &speciesName) const =0
 Specie mass diffusivity for specie in a multicomponent. More...
 
virtual tmp< volScalarFieldL (const word &speciesName, const volScalarField &Tf) const =0
 Latent heat (delta Hc) More...
 
virtual tmp< volScalarFieldKexp (const volScalarField &field)=0
 Explicit full mass transfer. More...
 
virtual tmp< volScalarFieldKSp (label modelVariable, const volScalarField &field)=0
 Implicit mass transfer. More...
 
virtual tmp< volScalarFieldKSu (label modelVariable, const volScalarField &field)=0
 Explicit mass transfer. More...
 
virtual const dimensionedScalarTactivate () const noexcept=0
 Reference value. More...
 
virtual bool includeDivU () const noexcept
 
bool includeVolChange ()
 Add volume change in pEq. More...
 
const wordvariable () 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< interfaceCompositionModelNew (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< volScalarFieldgetSpecieMassFraction (const word &speciesName, const pureMixture< ThermoType > &thermo) const
 Return mass fraction for a pureMixture equal to one. More...
 
template<class ThermoType >
tmp< volScalarFieldgetSpecieMassFraction (const word &speciesName, const multiComponentMixture< ThermoType > &thermo) const
 Return mass fraction for speciesName. More...
 
template<class ThermoType >
tmp< volScalarFieldMwMixture (const pureMixture< ThermoType > &thermo) const
 Return moleculas weight of the mixture for pureMixture [Kg/mol]. More...
 
template<class ThermoType >
tmp< volScalarFieldMwMixture (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 phasePairpair_
 Phase pair. More...
 
word speciesName_
 Names of the transferring specie. More...
 
const fvMeshmesh_
 Reference to mesh. More...
 
- Static Protected Attributes inherited from interfaceCompositionModel
static const Enum< modelVariablemodelVariableNames_
 Selection names for the modelVariable. More...
 

Detailed Description

template<class Thermo, class OtherThermo>
class Foam::meltingEvaporationModels::kineticGasEvaporation< Thermo, OtherThermo >

Considering the Hertz Knudsen formula, which gives the evaporation-condensation flux based on the kinetic theory for flat interface:

\[ Flux = C \sqrt{\frac{M}{2 \pi R T_{activate}}}(p - p_{sat}) \]

where:

\( Flux \) = mass flux rate [kg/s/m2]
\( M \) = molecular weight
\( T_{activate} \) = saturation temperature
\( C \) = accommodation coefficient
\( R \) = universal gas constant
\( p_{sat} \) = saturation pressure
\( p \) = vapor partial pressure

The Clapeyron-Clausius equation relates the pressure to the temperature for the saturation condition:

\[ \frac{dp}{dT} = - \frac{L}{T (\nu_v - \nu_l)} \]

where:

\( L \) = latent heat
\( \nu_v \) = inverse of the vapor density
\( \nu_l \) = inverse of the liquid density

Using the above relations:

\[ Flux = 2 \frac{C}{2 - C} \sqrt{\frac{M}{2 \pi R {T_activate}^3}} L (\rho_{v}) (T - T_{activate}) \]

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:

  • Van P. Carey, Liquid-Vapor Phase Change Phenomena, ISBN 0-89116836, 1992, pp. 112-121.
Usage

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
Source files

Definition at line 205 of file kineticGasEvaporation.H.

Constructor & Destructor Documentation

◆ kineticGasEvaporation()

◆ ~kineticGasEvaporation()

virtual ~kineticGasEvaporation ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "kineticGasEvaporation< Thermo, OtherThermo >"  )

Runtime type information.

◆ Kexp()

Foam::tmp< Foam::volScalarField > Kexp ( const volScalarField field)
virtual

Explicit total mass transfer coefficient.

Implements interfaceCompositionModel.

Definition at line 150 of file kineticGasEvaporation.C.

References Foam::dimDensity, Foam::dimTemperature, L, Foam::mag(), Foam::max(), IOobject::member(), mesh, Foam::constant::mathematical::pi(), Foam::pow3(), Foam::constant::physicoChemical::R, tmp< T >::ref(), Foam::sign(), Foam::sqrt(), T, T0, and Foam::Zero.

Here is the call graph for this function:

◆ KSp()

Foam::tmp< Foam::volScalarField > KSp ( label  modelVariable,
const volScalarField field 
)
virtual

Implicit mass transfer coefficient.

Implements interfaceCompositionModel.

Definition at line 225 of file kineticGasEvaporation.C.

References Foam::pos(), and Foam::sign().

Here is the call graph for this function:

◆ KSu()

Foam::tmp< Foam::volScalarField > KSu ( label  modelVariable,
const volScalarField field 
)
virtual

Explicit mass transfer coefficient.

Implements interfaceCompositionModel.

Definition at line 253 of file kineticGasEvaporation.C.

References Foam::pos(), and Foam::sign().

Here is the call graph for this function:

◆ Tactivate()

virtual const dimensionedScalar & Tactivate ( ) const
inlinevirtualnoexcept

Return Tactivate.

Implements interfaceCompositionModel.

Definition at line 282 of file kineticGasEvaporation.H.

◆ includeDivU()

virtual bool includeDivU ( ) const
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 289 of file kineticGasEvaporation.H.


The documentation for this class was generated from the following files: