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...

Collaboration diagram for interfaceCompositionModel:
[legend]

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...
 
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)
 

Protected Attributes

const phasePairpair_
 Phase pair. More...
 
const hashedWordList speciesNames_
 Names of the transferring species. 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 60 of file interfaceCompositionModel.H.

Constructor & Destructor Documentation

◆ interfaceCompositionModel()

interfaceCompositionModel ( const dictionary dict,
const phasePair pair 
)

Construct from a dictionary and a phase pair.

◆ ~interfaceCompositionModel()

virtual ~interfaceCompositionModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ TypeName()

TypeName ( "interfaceCompositionModel"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

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

◆ New()

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

◆ update()

virtual void update ( const volScalarField Tf)
pure virtual

Update the composition.

◆ species()

const hashedWordList & species ( ) const
inline

The transferring species names.

Definition at line 123 of file interfaceCompositionModel.H.

◆ 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()

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

Interface mass fraction.

◆ YfPrime()

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

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

◆ dY()

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

Mass fraction difference between the interface and the field.

◆ D()

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

Mass diffusivity.

◆ L()

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

Latent heat.

◆ addMDotL()

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

Add latent heat flow rate to total.

Member Data Documentation

◆ pair_

const phasePair& pair_
protected

Phase pair.

Definition at line 67 of file interfaceCompositionModel.H.

◆ speciesNames_

const hashedWordList speciesNames_
protected

Names of the transferring species.

Definition at line 70 of file interfaceCompositionModel.H.


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