Henry< Thermo, OtherThermo > Class Template Reference

Henry's law for gas solubility in liquid. The concentration of a dissolved species in the liquid is proportional to its partial pressure in the gas. A dimensionless solubility, \(k\), is given for each species. This is the ratio of the concentration of the species in the liquid to the corresponding concentration in the gas; i.e., \(k = c_{i,liq}/c_{i,gas}\). Mixing in the gas is assumed to be ideal. More...

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

Public Member Functions

 TypeName ("Henry")
 Runtime type information. More...
 
 Henry (const dictionary &dict, const phasePair &pair)
 Construct from components. More...
 
virtual ~Henry ()=default
 Destructor. More...
 
virtual void update (const volScalarField &Tf)
 Update the composition. More...
 
virtual tmp< volScalarFieldYf (const word &speciesName, const volScalarField &Tf) const
 The interface species fraction. More...
 
virtual tmp< volScalarFieldYfPrime (const word &speciesName, const volScalarField &Tf) const
 The interface species fraction derivative w.r.t. temperature. 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< 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::interfaceCompositionModels::Henry< Thermo, OtherThermo >

Henry's law for gas solubility in liquid. The concentration of a dissolved species in the liquid is proportional to its partial pressure in the gas. A dimensionless solubility, \(k\), is given for each species. This is the ratio of the concentration of the species in the liquid to the corresponding concentration in the gas; i.e., \(k = c_{i,liq}/c_{i,gas}\). Mixing in the gas is assumed to be ideal.

Source files

Definition at line 61 of file Henry.H.

Constructor & Destructor Documentation

◆ Henry()

Henry ( const dictionary dict,
const phasePair pair 
)

Construct from components.

Definition at line 33 of file Henry.C.

References Foam::dimless, Foam::exit(), Foam::FatalError, FatalErrorInFunction, Foam::name(), and UList< T >::size().

Here is the call graph for this function:

◆ ~Henry()

virtual ~Henry ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "Henry< Thermo, OtherThermo >"  )

Runtime type information.

◆ update()

void update ( const volScalarField Tf)
virtual

Update the composition.

Definition at line 65 of file Henry.C.

◆ Yf()

Foam::tmp< Foam::volScalarField > Yf ( const word speciesName,
const volScalarField Tf 
) const
virtual

The interface species fraction.

Reimplemented from InterfaceCompositionModel< Thermo, OtherThermo >.

Definition at line 81 of file Henry.C.

◆ YfPrime()

Foam::tmp< Foam::volScalarField > YfPrime ( const word speciesName,
const volScalarField Tf 
) const
virtual

The interface species fraction derivative w.r.t. temperature.

Definition at line 108 of file Henry.C.

References Foam::dimless, Foam::dimTemperature, IOobject::groupName(), and Time::New().

Here is the call graph for this function:

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