interfaceCompositionModel Class Referenceabstract

Generic base class for interface models. Mass transfer models are interface models between two thermos. Abstract class for mass transfer functions. More...

Inheritance diagram for interfaceCompositionModel:
[legend]
Collaboration diagram for interfaceCompositionModel:
[legend]

Public Types

enum  modelVariable { T, P, Y, alpha }
 Enumeration for variable based mass transfer models. More...
 

Public Member Functions

 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...
 
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< volScalarFieldD (const word &speciesName) const =0
 Mass diffusivity. 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...
 
 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...
 
virtual void update (const volScalarField &Tf)=0
 Update the composition. More...
 
const hashedWordListspecies () const
 The transferring species names. More...
 
bool transports (word &speciesName) const
 
virtual tmp< volScalarFieldYf (const word &speciesName, const volScalarField &Tf) const =0
 Interface mass fraction. More...
 
virtual tmp< volScalarFieldYfPrime (const word &speciesName, const volScalarField &Tf) const =0
 The interface mass fraction derivative w.r.t. temperature. More...
 
virtual tmp< volScalarFielddY (const word &speciesName, const volScalarField &Tf) const =0
 Mass fraction difference between the interface and the field. More...
 
virtual tmp< volScalarFieldD (const word &speciesName) const =0
 Mass diffusivity. More...
 
virtual tmp< volScalarFieldL (const word &speciesName, const volScalarField &Tf) const =0
 Latent heat. More...
 
virtual void addMDotL (const volScalarField &K, const volScalarField &Tf, volScalarField &mDotL, volScalarField &mDotLPrime) const =0
 Add latent heat flow rate to total. More...
 

Static Public Member Functions

static autoPtr< interfaceCompositionModelNew (const dictionary &dict, const phasePair &pair)
 
static autoPtr< interfaceCompositionModelNew (const dictionary &dict, const phasePair &pair)
 

Protected Attributes

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...
 
const hashedWordList speciesNames_
 Names of the transferring species. More...
 

Static Protected Attributes

static const Enum< modelVariablemodelVariableNames_
 Selection names for the modelVariable. More...
 

Detailed Description

Generic base class for interface models. Mass transfer models are interface models between two thermos. Abstract class for mass transfer functions.

Generic base class for interface composition models. These models describe the composition in phase 1 of the supplied pair at the interface with phase 2.

Source files

Source files

Definition at line 59 of file interfaceCompositionModel.H.

Member Enumeration Documentation

◆ modelVariable

Enumeration for variable based mass transfer models.

Enumerator
alpha 

Definition at line 66 of file interfaceCompositionModel.H.

Constructor & Destructor Documentation

◆ interfaceCompositionModel() [1/2]

interfaceCompositionModel ( const dictionary dict,
const phasePair pair 
)

Construct from a dictionary and a phase pair.

Definition at line 54 of file interfaceCompositionModel.C.

◆ ~interfaceCompositionModel() [1/2]

virtual ~interfaceCompositionModel ( )
virtualdefault

Destructor.

◆ interfaceCompositionModel() [2/2]

interfaceCompositionModel ( const dictionary dict,
const phasePair pair 
)

Construct from a dictionary and a phase pair.

◆ ~interfaceCompositionModel() [2/2]

virtual ~interfaceCompositionModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName() [1/2]

TypeName ( "interfaceCompositionModel"  )

Runtime type information.

◆ declareRunTimeSelectionTable() [1/2]

declareRunTimeSelectionTable ( autoPtr  ,
interfaceCompositionModel  ,
dictionary  ,
(const dictionary &dict, const phasePair &pair ,
(dict, pair  
)

◆ New() [1/2]

Foam::autoPtr< Foam::interfaceCompositionModel > New ( const dictionary dict,
const phasePair pair 
)
static

Definition at line 79 of file interfaceCompositionModel.C.

References dict, Foam::endl(), Foam::exit(), Foam::FatalIOError, FatalIOErrorInLookup, Foam::Info, phasePair::phase1(), phasePair::phase2(), and phaseModel::thermo().

Here is the call graph for this function:

◆ transferSpecie()

const Foam::word transferSpecie ( ) const

Return the transferring species name.

Definition at line 116 of file interfaceCompositionModel.C.

References interfaceCompositionModel::speciesName_.

Referenced by MassTransferPhaseSystem< BasePhaseSystem >::calculateL().

Here is the caller graph for this function:

◆ pair()

const Foam::phasePair & pair ( ) const

The phase pair.

Definition at line 122 of file interfaceCompositionModel.C.

◆ Yf() [1/2]

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

Interface mass fraction.

Implemented in InterfaceCompositionModel< Thermo, OtherThermo >, NonRandomTwoLiquid< Thermo, OtherThermo >, Saturated< Thermo, OtherThermo >, Henry< Thermo, OtherThermo >, and Raoult< Thermo, OtherThermo >.

Referenced by InterfaceCompositionPhaseChangePhaseSystem< BasePhaseSystem >::massTransfer().

Here is the caller graph for this function:

◆ dY() [1/2]

virtual tmp<volScalarField> dY ( const word speciesName,
const volScalarField Tf 
) const
pure virtual

Mass fraction difference between the interface and the field.

Implemented in InterfaceCompositionModel< Thermo, OtherThermo >, and InterfaceCompositionModel< Thermo, OtherThermo >.

◆ D() [1/2]

virtual tmp<volScalarField> D ( const word speciesName) const
pure virtual

Mass diffusivity.

Implemented in InterfaceCompositionModel< Thermo, OtherThermo >, and InterfaceCompositionModel< Thermo, OtherThermo >.

Referenced by InterfaceCompositionPhaseChangePhaseSystem< BasePhaseSystem >::massTransfer().

Here is the caller graph for this function:

◆ L() [1/2]

virtual tmp<volScalarField> L ( const word speciesName,
const volScalarField Tf 
) const
pure virtual

Latent heat (delta Hc)

Implemented in InterfaceCompositionModel< Thermo, OtherThermo >, and InterfaceCompositionModel< Thermo, OtherThermo >.

Referenced by MassTransferPhaseSystem< BasePhaseSystem >::calculateL().

Here is the caller graph for this function:

◆ Kexp()

virtual tmp<volScalarField> Kexp ( const volScalarField field)
pure virtual

Explicit full mass transfer.

Implemented in interfaceOxideRate< Thermo, OtherThermo >, kineticGasEvaporation< Thermo, OtherThermo >, interfaceHeatResistance< Thermo, OtherThermo >, and Lee< Thermo, OtherThermo >.

Referenced by MassTransferPhaseSystem< BasePhaseSystem >::correctMassSources().

Here is the caller graph for this function:

◆ KSp()

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

Implicit mass transfer.

Implemented in interfaceOxideRate< Thermo, OtherThermo >, kineticGasEvaporation< Thermo, OtherThermo >, interfaceHeatResistance< Thermo, OtherThermo >, and Lee< Thermo, OtherThermo >.

Referenced by MassTransferPhaseSystem< BasePhaseSystem >::heatTransfer(), and MassTransferPhaseSystem< BasePhaseSystem >::volTransfer().

Here is the caller graph for this function:

◆ KSu()

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

◆ Tactivate()

◆ includeDivU()

bool includeDivU ( ) const
virtualnoexcept

Add/subtract alpha*div(U) as a source term for alpha, substituting div(U) = mDot(1/rho1 - 1/rho2)

Reimplemented in interfaceOxideRate< Thermo, OtherThermo >, kineticGasEvaporation< Thermo, OtherThermo >, interfaceHeatResistance< Thermo, OtherThermo >, and Lee< Thermo, OtherThermo >.

Definition at line 134 of file interfaceCompositionModel.C.

Referenced by MassTransferPhaseSystem< BasePhaseSystem >::alphaTransfer().

Here is the caller graph for this function:

◆ includeVolChange()

bool includeVolChange ( )

Add volume change in pEq.

Definition at line 140 of file interfaceCompositionModel.C.

◆ variable()

const Foam::word & variable ( ) const

Returns the variable on which the model is based.

Definition at line 128 of file interfaceCompositionModel.C.

◆ TypeName() [2/2]

TypeName ( "interfaceCompositionModel"  )

Runtime type information.

◆ declareRunTimeSelectionTable() [2/2]

declareRunTimeSelectionTable ( autoPtr  ,
interfaceCompositionModel  ,
dictionary  ,
(const dictionary &dict, const phasePair &pair ,
(dict, pair  
)

◆ New() [2/2]

static autoPtr<interfaceCompositionModel> New ( const dictionary dict,
const phasePair pair 
)
static

◆ update()

virtual void update ( const volScalarField Tf)
pure virtual

◆ species()

const hashedWordList& species ( ) const
inline

The transferring species names.

Definition at line 123 of file interfaceCompositionModel.H.

References interfaceCompositionModel::speciesNames_.

Referenced by InterfaceCompositionPhaseChangePhaseSystem< BasePhaseSystem >::dmdts(), InterfaceCompositionPhaseChangePhaseSystem< BasePhaseSystem >::InterfaceCompositionPhaseChangePhaseSystem(), and InterfaceCompositionPhaseChangePhaseSystem< BasePhaseSystem >::massTransfer().

Here is the caller graph for this function:

◆ transports()

bool transports ( word speciesName) const

Returns whether the species is transported by the model and provides the name of the diffused species

Definition at line 96 of file interfaceCompositionModel.C.

◆ Yf() [2/2]

◆ YfPrime()

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

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

Implemented in NonRandomTwoLiquid< Thermo, OtherThermo >, Saturated< Thermo, OtherThermo >, Henry< Thermo, OtherThermo >, and Raoult< Thermo, OtherThermo >.

◆ dY() [2/2]

virtual tmp<volScalarField> dY ( const word speciesName,
const volScalarField Tf 
) const
pure virtual

Mass fraction difference between the interface and the field.

Implemented in InterfaceCompositionModel< Thermo, OtherThermo >, and InterfaceCompositionModel< Thermo, OtherThermo >.

◆ D() [2/2]

virtual tmp<volScalarField> D ( const word speciesName) const
pure virtual

◆ L() [2/2]

virtual tmp<volScalarField> L ( const word speciesName,
const volScalarField Tf 
) const
pure virtual

◆ addMDotL()

virtual void addMDotL ( const volScalarField K,
const volScalarField Tf,
volScalarField mDotL,
volScalarField mDotLPrime 
) const
pure virtual

Add latent heat flow rate to total.

Implemented in InterfaceCompositionModel< Thermo, OtherThermo >.

Member Data Documentation

◆ modelVariableNames_

const Foam::Enum< Foam::interfaceCompositionModel::modelVariable > modelVariableNames_
staticprotected
Initial value:
{
{ modelVariable::T, "temperature" },
{ modelVariable::P, "pressure" },
{ modelVariable::Y, "massFraction" },
{ modelVariable::alpha, "alphaVolumeFraction" },
}

Selection names for the modelVariable.

Definition at line 79 of file interfaceCompositionModel.H.

◆ modelVariable_

modelVariable modelVariable_
protected

Enumeration for the model variable.

Definition at line 82 of file interfaceCompositionModel.H.

◆ includeVolChange_

bool includeVolChange_
protected

Add volume change in pEq.

Definition at line 85 of file interfaceCompositionModel.H.

◆ pair_

const phasePair & pair_
protected

Phase pair.

Definition at line 88 of file interfaceCompositionModel.H.

◆ speciesName_

word speciesName_
protected

Names of the transferring specie.

Definition at line 91 of file interfaceCompositionModel.H.

Referenced by interfaceCompositionModel::transferSpecie().

◆ mesh_

const fvMesh& mesh_
protected

Reference to mesh.

Definition at line 94 of file interfaceCompositionModel.H.

◆ speciesNames_

const hashedWordList speciesNames_
protected

Names of the transferring species.

Definition at line 70 of file interfaceCompositionModel.H.

Referenced by interfaceCompositionModel::species().


The documentation for this class was generated from the following files:
Foam::constant::atomic::alpha
const dimensionedScalar alpha
Fine-structure constant: default SI units: [].
Definition: readThermalProperties.H:212
T
const volScalarField & T
Definition: createFieldRefs.H:2
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7