PhaseChangeModel< CloudType > Class Template Referenceabstract

Templated phase change model class. More...

Inheritance diagram for PhaseChangeModel< CloudType >:
[legend]
Collaboration diagram for PhaseChangeModel< CloudType >:
[legend]

Public Types

enum  enthalpyTransferType { etLatentHeat, etEnthalpyDifference }
 Enthalpy transfer type. More...
 

Public Member Functions

 TypeName ("phaseChangeModel")
 Runtime type information. More...
 
 declareRunTimeSelectionTable (autoPtr, PhaseChangeModel, dictionary,(const dictionary &dict, CloudType &owner),(dict, owner))
 Declare runtime constructor selection table. More...
 
 PhaseChangeModel (CloudType &owner)
 Construct null from owner. More...
 
 PhaseChangeModel (const dictionary &dict, CloudType &owner, const word &type)
 Construct from dictionary. More...
 
 PhaseChangeModel (const PhaseChangeModel< CloudType > &pcm)
 Construct copy. More...
 
virtual autoPtr< PhaseChangeModel< CloudType > > clone () const =0
 Construct and return a clone. More...
 
virtual ~PhaseChangeModel ()=default
 Destructor. More...
 
const enthalpyTransferTypeenthalpyTransfer () const
 Return the enthalpy transfer type enumeration. More...
 
virtual void calculate (const scalar dt, const label celli, const scalar Re, const scalar Pr, const scalar d, const scalar nu, const scalar rho, const scalar T, const scalar Ts, const scalar pc, const scalar Tc, const scalarField &X, const scalarField &solMass, const scalarField &liqMass, scalarField &dMassPC) const =0
 Update model. More...
 
virtual scalar dh (const label idc, const label idl, const scalar p, const scalar T) const
 Return the enthalpy per unit mass. More...
 
virtual scalar Tvap (const scalarField &X) const
 Return vapourisation temperature. More...
 
virtual scalar TMax (const scalar p, const scalarField &X) const
 Return maximum/limiting temperature. More...
 
void addToPhaseChangeMass (const scalar dMass)
 Add to phase change mass. More...
 
virtual void info (Ostream &os)
 Write injection info to stream. More...
 

Static Public Member Functions

static autoPtr< PhaseChangeModel< CloudType > > New (const dictionary &dict, CloudType &owner)
 Selector. More...
 

Static Public Attributes

static const wordList enthalpyTransferTypeNames
 Name representations of enthalpy transfer types. More...
 

Protected Member Functions

enthalpyTransferType wordToEnthalpyTransfer (const word &etName) const
 Convert word to enthalpy transfer type. More...
 
scalar Sh () const
 Sherwood number. More...
 

Protected Attributes

enthalpyTransferType enthalpyTransfer_
 Enthalpy transfer type enumeration. More...
 
scalar dMass_
 Mass of lagrangian phase converted. More...
 

Detailed Description

template<class CloudType>
class Foam::PhaseChangeModel< CloudType >

Templated phase change model class.

Source files

Definition at line 61 of file ReactingCloud.H.

Member Enumeration Documentation

◆ enthalpyTransferType

Enthalpy transfer type.

Enumerator
etLatentHeat 
etEnthalpyDifference 

Definition at line 68 of file PhaseChangeModel.H.

Constructor & Destructor Documentation

◆ PhaseChangeModel() [1/3]

Construct null from owner.

Definition at line 67 of file PhaseChangeModel.C.

◆ PhaseChangeModel() [2/3]

PhaseChangeModel ( const dictionary dict,
CloudType owner,
const word type 
)

Construct from dictionary.

Definition at line 91 of file PhaseChangeModel.C.

◆ PhaseChangeModel() [3/3]

Construct copy.

Definition at line 79 of file PhaseChangeModel.C.

◆ ~PhaseChangeModel()

virtual ~PhaseChangeModel ( )
virtualdefault

Destructor.

Member Function Documentation

◆ wordToEnthalpyTransfer()

Foam::PhaseChangeModel< CloudType >::enthalpyTransferType wordToEnthalpyTransfer ( const word etName) const
protected

Convert word to enthalpy transfer type.

Definition at line 44 of file PhaseChangeModel.C.

◆ Sh()

scalar Sh ( ) const
protected

Sherwood number.

◆ TypeName()

TypeName ( "phaseChangeModel"  )

Runtime type information.

◆ declareRunTimeSelectionTable()

declareRunTimeSelectionTable ( autoPtr  ,
PhaseChangeModel< CloudType ,
dictionary  ,
(const dictionary &dict, CloudType &owner)  ,
(dict, owner)   
)

Declare runtime constructor selection table.

◆ clone()

◆ New()

Foam::autoPtr< Foam::PhaseChangeModel< CloudType > > New ( const dictionary dict,
CloudType owner 
)
static

Selector.

Definition at line 36 of file PhaseChangeModelNew.C.

◆ enthalpyTransfer()

const Foam::PhaseChangeModel< CloudType >::enthalpyTransferType & enthalpyTransfer ( ) const

Return the enthalpy transfer type enumeration.

Definition at line 110 of file PhaseChangeModel.C.

◆ calculate()

virtual void calculate ( const scalar  dt,
const label  celli,
const scalar  Re,
const scalar  Pr,
const scalar  d,
const scalar  nu,
const scalar  rho,
const scalar  T,
const scalar  Ts,
const scalar  pc,
const scalar  Tc,
const scalarField X,
const scalarField solMass,
const scalarField liqMass,
scalarField dMassPC 
) const
pure virtual

◆ dh()

Foam::scalar dh ( const label  idc,
const label  idl,
const scalar  p,
const scalar  T 
) const
virtual

Return the enthalpy per unit mass.

Reimplemented in LiquidEvapFuchsKnudsen< CloudType >, LiquidEvaporationBoil< CloudType >, and LiquidEvaporation< CloudType >.

Definition at line 118 of file PhaseChangeModel.C.

◆ Tvap()

Foam::scalar Tvap ( const scalarField X) const
virtual

Return vapourisation temperature.

Reimplemented in LiquidEvapFuchsKnudsen< CloudType >, LiquidEvaporationBoil< CloudType >, and LiquidEvaporation< CloudType >.

Definition at line 141 of file PhaseChangeModel.C.

◆ TMax()

Foam::scalar TMax ( const scalar  p,
const scalarField X 
) const
virtual

Return maximum/limiting temperature.

Reimplemented in LiquidEvapFuchsKnudsen< CloudType >, LiquidEvaporationBoil< CloudType >, and LiquidEvaporation< CloudType >.

Definition at line 131 of file PhaseChangeModel.C.

◆ addToPhaseChangeMass()

void addToPhaseChangeMass ( const scalar  dMass)

Add to phase change mass.

Definition at line 148 of file PhaseChangeModel.C.

◆ info()

void info ( Ostream os)
virtual

Write injection info to stream.

Definition at line 155 of file PhaseChangeModel.C.

Member Data Documentation

◆ enthalpyTransferTypeNames

const Foam::wordList enthalpyTransferTypeNames
static
Initial value:
{
"latentHeat", "enthalpyDifference"
}

Name representations of enthalpy transfer types.

Definition at line 75 of file PhaseChangeModel.H.

◆ enthalpyTransfer_

enthalpyTransferType enthalpyTransfer_
protected

Enthalpy transfer type enumeration.

Definition at line 83 of file PhaseChangeModel.H.

◆ dMass_

scalar dMass_
protected

Mass of lagrangian phase converted.

Definition at line 89 of file PhaseChangeModel.H.


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