Mass transfer Lee model. Simple model driven by field value difference as: More...
#include <Lee.H>


Public Member Functions | |
| TypeName ("Lee") | |
| Runtime type information. | |
| Lee (const dictionary &dict, const phasePair &pair) | |
| Construct from components. | |
| virtual | ~Lee ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | Kexp (const volScalarField &field) |
| Explicit total mass transfer coefficient. | |
| virtual tmp< volScalarField > | KSp (label modelVariable, const volScalarField &field) |
| Implicit mass transfer coefficient. | |
| virtual tmp< volScalarField > | KSu (label modelVariable, const volScalarField &field) |
| Explicit mass transfer coefficient. | |
| virtual const dimensionedScalar & | Tactivate () const noexcept |
| Return T transition between phases. | |
| virtual bool | includeDivU () const noexcept |
| Add/subtract alpha*div(U) as a source term for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2). | |
| Public Member Functions inherited from InterfaceCompositionModel< Thermo, OtherThermo > | |
| InterfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
| Construct from components. | |
| virtual | ~InterfaceCompositionModel ()=default |
| Destructor. | |
| virtual tmp< volScalarField > | dY (const word &speciesName, const volScalarField &Tf) const |
| Mass fraction difference between the interface and the field. | |
| virtual tmp< volScalarField > | Yf (const word &speciesName, const volScalarField &Tf) const |
| Reference mass fraction for species based models. | |
| virtual tmp< volScalarField > | Dfrom (const word &speciesName) const |
| Specie mass diffusivity for pure mixture. | |
| virtual tmp< volScalarField > | Dto (const word &speciesName) const |
| Specie mass diffusivity for specie in a multicomponent. | |
| virtual tmp< volScalarField > | L (const word &speciesName, const volScalarField &Tf) const |
| Latent heat (to - from)(thermo - otherThermo). | |
| InterfaceCompositionModel (const dictionary &dict, const phasePair &pair) | |
| Construct from components. | |
| ~InterfaceCompositionModel ()=default | |
| Destructor. | |
| virtual tmp< volScalarField > | dY (const word &speciesName, const volScalarField &Tf) const |
| Mass fraction difference between the interface and the field. | |
| virtual tmp< volScalarField > | D (const word &speciesName) const |
| Mass diffusivity. | |
| virtual tmp< volScalarField > | L (const word &speciesName, const volScalarField &Tf) const |
| Latent heat. | |
| virtual void | addMDotL (const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const |
| Add latent heat flow rate to total. | |
| 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. | |
| 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. | |
| virtual | ~interfaceCompositionModel ()=default |
| Destructor. | |
| const word | transferSpecie () const |
| Return the transferring species name. | |
| const phasePair & | pair () const |
| The phase pair. | |
| const multiphaseInterSystem & | fluid () const |
| Return the multiphaseInterSystem this interface belongs to. | |
| bool | includeVolChange () |
| Add volume change in pEq. | |
| const word & | variable () const |
| Returns the variable on which the model is based. | |
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. | |
| 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. | |
| template<class ThermoType> | |
| tmp< volScalarField > | getSpecieMassFraction (const word &speciesName, const pureMixture< ThermoType > &thermo) const |
| Return mass fraction for a pureMixture equal to one. | |
| template<class ThermoType> | |
| tmp< volScalarField > | getSpecieMassFraction (const word &speciesName, const multiComponentMixture< ThermoType > &thermo) const |
| Return mass fraction for speciesName. | |
| template<class ThermoType> | |
| tmp< volScalarField > | MwMixture (const pureMixture< ThermoType > &thermo) const |
| Return moleculas weight of the mixture for pureMixture [Kg/mol]. | |
| template<class ThermoType> | |
| tmp< volScalarField > | MwMixture (const multiComponentMixture< ThermoType > &) const |
| Return moleculas weight of the mixture for multiComponentMixture. | |
| 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. | |
| 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. | |
| Protected Attributes inherited from InterfaceCompositionModel< Thermo, OtherThermo > | |
| const Thermo & | fromThermo_ |
| Thermo (from). | |
| const OtherThermo & | toThermo_ |
| Other Thermo (to). | |
| const dimensionedScalar | Le_ |
| Lewis number. | |
| const Thermo & | thermo_ |
| Thermo. | |
| const OtherThermo & | otherThermo_ |
| Other Thermo. | |
| Protected Attributes inherited from interfaceCompositionModel | |
| modelVariable | modelVariable_ |
| Enumeration for the model variable. | |
| bool | includeVolChange_ |
| Add volume change in pEq. | |
| const phasePair & | pair_ |
| Phase pair. | |
| word | speciesName_ |
| Names of the transferring specie. | |
| const fvMesh & | mesh_ |
| Reference to mesh. | |
| Static Protected Attributes inherited from interfaceCompositionModel | |
| static const Enum< modelVariable > | modelVariableNames_ |
| Selection names for the modelVariable. | |
Mass transfer Lee model. Simple model driven by field value difference as:
![\[ \dot{m} = C \rho \alpha (T - T_{activate})/T_{activate}
\]](form_763.png)
where C is a model constant.
if C > 0:
![\[ \dot{m} = C \rho \alpha (T - T_{activate})/T_{activate}
\]](form_763.png)
for
![\[ T > T_{activate} \]](form_764.png)
and
![\[ mDot = 0.0 \]](form_765.png)
for
![\[ T < T_{activate} \]](form_766.png)
if C < 0:
![\[ \dot{m} = -C \rho \alpha (T_{activate} - T)/T_{activate}
\]](form_767.png)
for
![\[ T < T_{activate} \]](form_768.png)
and
![\[ \dot{m} = 0.0 \]](form_769.png)
for
![\[ T > T_{activate} \]](form_764.png)
Based on the reference:
massTransferModel
(
(solid to liquid)
{
type Lee;
C 40;
Tactivate 302.78;
}
);
Where:
| Property | Description | Required | Default value |
|---|---|---|---|
Tactivate | Activation temperature | yes | |
C | Model constant | yes | |
includeVolChange | Volumen change | no | yes |
species | Specie name on the other phase | no | none |
| Lee | ( | const dictionary & | dict, |
| const phasePair & | pair ) |
Construct from components.
Definition at line 28 of file Lee.C.
References dict, Foam::dimTemperature, Foam::dimTime, InterfaceCompositionModel< Thermo, OtherThermo >::InterfaceCompositionModel(), Foam::inv(), and interfaceCompositionModel::pair().

|
virtualdefault |
| TypeName | ( | "Lee< Thermo, OtherThermo >" | ) |
Runtime type information.
References dict, and interfaceCompositionModel::pair().

|
virtual |
Explicit total mass transfer coefficient.
Implements interfaceCompositionModel.
Definition at line 45 of file Lee.C.
References Foam::clamp(), interfaceCompositionModel::pair(), Foam::pos(), rho, and Foam::sign().

|
virtual |
Implicit mass transfer coefficient.
Implements interfaceCompositionModel.
Definition at line 83 of file Lee.C.
References Foam::clamp(), interfaceCompositionModel::modelVariable_, interfaceCompositionModel::pair(), Foam::pos(), rho, Foam::sign(), and interfaceCompositionModel::variable().

|
virtual |
Explicit mass transfer coefficient.
Implements interfaceCompositionModel.
Definition at line 126 of file Lee.C.
References Foam::clamp(), interfaceCompositionModel::modelVariable_, interfaceCompositionModel::pair(), Foam::pos(), rho, Foam::sign(), and interfaceCompositionModel::variable().

|
inlinevirtualnoexcept |
Return T transition between phases.
Implements interfaceCompositionModel.
Definition at line 206 of file Lee.H.
References Foam::noexcept.
|
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 215 of file Lee.H.
References Foam::noexcept.