FreeFOAM The Cross-Platform CFD Toolkit
Hosted by SourceForge:
Get FreeFOAM at SourceForge.net.
            Fast, secure and Free Open Source software downloads

CompositionModel< CloudType > Class Template Reference

Templated reacting parcel composition model class Consists of carrier species (via thermo package), and additional liquids and solids. More...

#include <lagrangianIntermediate/CompositionModel.H>


Detailed Description

template<class CloudType>
class Foam::CompositionModel< CloudType >

Templated reacting parcel composition model class Consists of carrier species (via thermo package), and additional liquids and solids.

Source files

Definition at line 63 of file CompositionModel.H.

Inheritance diagram for CompositionModel< CloudType >:
Collaboration diagram for CompositionModel< CloudType >:

List of all members.

Public Member Functions

 TypeName ("CompositionModel")
 Runtime type information.
 declareRunTimeSelectionTable (autoPtr, CompositionModel, dictionary,(const dictionary &dict, CloudType &owner),(dict, owner))
 Declare runtime constructor selection table.
 CompositionModel (const dictionary &dict, CloudType &owner, const word &type)
 Construct from dictionary.
virtual  ~CompositionModel ()
 Destructor.
const CloudType &  owner () const
 Return the cloud object.
const dictionary &  dict () const
 Return the cloud dictionary.
const dictionary &  coeffDict () const
 Return the coefficients dictionary.
const multiComponentMixture
< typename
CloudType::thermoType > &  
mcCarrierThermo () const
 Return the carrier phase thermo package.
const liquidMixture &  liquids () const
 Return the global (additional) liquids.
const solidMixture &  solids () const
 Return the global (additional) solids.
const phasePropertiesList &  phaseProps () const
 Return the list of phase properties.
label  nPhase () const
 Return the number of phases.
const wordList &  phaseTypes () const
 Return the list of phase type names.
const wordList &  stateLabels () const
 Return the list of state labels (s), (l), (g) etc.
const wordList &  componentNames (const label phaseI) const
 Return the list of component names for phaseI.
label  globalCarrierId (const word &cmptName) const
 Return global id of component cmptName in carrier thermo.
label  globalId (const label phaseI, const word &cmptName) const
 Return global id of component cmptName in phase phaseI.
const labelList &  globalIds (const label phaseI) const
 Return global ids of for phase phaseI.
label  localId (const label phaseI, const word &cmptName) const
 Return local id of component cmptName in phase phaseI.
label  localToGlobalCarrierId (const label phaseI, const label id) const
 Return global carrier id of component given local id.
const scalarField &  Y0 (const label phaseI) const
 Return the list of phase phaseI mass fractions.
scalarField  X (const label phaseI, const scalarField &Y) const
 Return the list of phase phaseI volume fractions fractions.
virtual const scalarField &  YMixture0 () const =0
 Return the list of mixture mass fractions.
virtual label  idGas () const =0
 Gas id.
virtual label  idLiquid () const =0
 Liquid id.
virtual label  idSolid () const =0
 Solid id.
virtual scalar  H (const label phaseI, const scalarField &Y, const scalar p, const scalar T) const
 Return total enthalpy for the phase phaseI.
virtual scalar  Hs (const label phaseI, const scalarField &Y, const scalar p, const scalar T) const
 Return sensible enthalpy for the phase phaseI.
virtual scalar  Hc (const label phaseI, const scalarField &Y, const scalar p, const scalar T) const
 Return chemical enthalpy for the phase phaseI.
virtual scalar  cp (const label phaseI, const scalarField &Y, const scalar p, const scalar T) const
 Return specific heat caoacity for the phase phaseI.
virtual scalar  L (const label phaseI, const scalarField &Y, const scalar p, const scalar T) const
 Return latent heat for the phase phaseI.

Static Public Member Functions

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

Constructor & Destructor Documentation

CompositionModel ( const dictionary &   dict,
CloudType &   owner,
const word &   type  
)

Construct from dictionary.

Definition at line 32 of file CompositionModel.C.

~CompositionModel (  ) [virtual]

Destructor.

Definition at line 75 of file CompositionModel.C.


Member Function Documentation

TypeName ( "CompositionModel< CloudType >"    )

Runtime type information.

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

Declare runtime constructor selection table.

Foam::autoPtr< Foam::CompositionModel< CloudType > > New ( const dictionary &   dict,
CloudType &   owner  
) [static]
const CloudType & owner (  ) const

Return the cloud object.

Definition at line 82 of file CompositionModel.C.

const Foam::dictionary & dict (  ) const

Return the cloud dictionary.

Definition at line 89 of file CompositionModel.C.

const Foam::dictionary & coeffDict (  ) const

Return the coefficients dictionary.

Definition at line 96 of file CompositionModel.C.

const Foam::multiComponentMixture< typename CloudType::thermoType > & mcCarrierThermo (  ) const

Return the carrier phase thermo package.

Definition at line 104 of file CompositionModel.C.

const Foam::liquidMixture & liquids (  ) const

Return the global (additional) liquids.

Definition at line 111 of file CompositionModel.C.

const Foam::solidMixture & solids (  ) const

Return the global (additional) solids.

Definition at line 118 of file CompositionModel.C.

const Foam::phasePropertiesList & phaseProps (  ) const

Return the list of phase properties.

Definition at line 126 of file CompositionModel.C.

Foam::label nPhase (  ) const

Return the number of phases.

Definition at line 133 of file CompositionModel.C.

Referenced by ReactingParcel< ParcelType >::readFields(), and ReactingParcel< ParcelType >::writeFields().

const Foam::wordList & phaseTypes (  ) const

Return the list of phase type names.

If only 1 phase, return the component names of that phase

Definition at line 141 of file CompositionModel.C.

Referenced by ReactingParcel< ParcelType >::ReactingParcel(), ReactingParcel< ParcelType >::readFields(), and ReactingParcel< ParcelType >::writeFields().

const Foam::wordList & componentNames ( const label   phaseI  ) const
Foam::label globalCarrierId ( const word &   cmptName  ) const

Return global id of component cmptName in carrier thermo.

Definition at line 173 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, forAll, and Foam::nl.

Foam::label globalId ( const label   phaseI,
const word &   cmptName  
) const

Return global id of component cmptName in phase phaseI.

Definition at line 200 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, and Foam::nl.

const Foam::labelList & globalIds ( const label   phaseI  ) const

Return global ids of for phase phaseI.

Definition at line 226 of file CompositionModel.C.

Foam::label localId ( const label   phaseI,
const word &   cmptName  
) const

Return local id of component cmptName in phase phaseI.

Definition at line 236 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, and Foam::nl.

Foam::label localToGlobalCarrierId ( const label   phaseI,
const label   id  
) const

Return global carrier id of component given local id.

Definition at line 262 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, and Foam::nl.

const Foam::scalarField & Y0 ( const label   phaseI  ) const

Return the list of phase phaseI mass fractions.

Definition at line 290 of file CompositionModel.C.

Foam::scalarField X ( const label   phaseI,
const scalarField &   Y  
) const

Return the list of phase phaseI volume fractions fractions.

based on supplied mass fractions Y

Definition at line 300 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, forAll, phaseProperties::globalIds(), Foam::nl, phaseProperties::phase(), and List< T >::size().

virtual const scalarField& YMixture0 (  ) const [pure virtual]

Return the list of mixture mass fractions.

If only 1 phase, return component fractions of that phase

Implemented in SingleMixtureFraction< CloudType >, and SinglePhaseMixture< CloudType >.

Foam::scalar H ( const label   phaseI,
const scalarField &   Y,
const scalar   p,
const scalar   T  
) const [virtual]

Return total enthalpy for the phase phaseI.

Definition at line 350 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, forAll, phaseProperties::globalIds(), Foam::nl, phaseProperties::phase(), and T.

Foam::scalar Hs ( const label   phaseI,
const scalarField &   Y,
const scalar   p,
const scalar   T  
) const [virtual]

Return sensible enthalpy for the phase phaseI.

Definition at line 414 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, forAll, phaseProperties::globalIds(), Foam::nl, phaseProperties::phase(), and T.

Foam::scalar Hc ( const label   phaseI,
const scalarField &   Y,
const scalar   p,
const scalar   T  
) const [virtual]

Return chemical enthalpy for the phase phaseI.

Definition at line 478 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, forAll, phaseProperties::globalIds(), Foam::nl, and phaseProperties::phase().

Foam::scalar cp ( const label   phaseI,
const scalarField &   Y,
const scalar   p,
const scalar   T  
) const [virtual]

Return specific heat caoacity for the phase phaseI.

Definition at line 538 of file CompositionModel.C.

References Foam::abort(), Foam::FatalError, FatalErrorIn, forAll, phaseProperties::globalIds(), Foam::nl, and phaseProperties::phase().

Foam::scalar L ( const label   phaseI,
const scalarField &   Y,
const scalar   p,
const scalar   T  
) const [virtual]

Return latent heat for the phase phaseI.

Definition at line 597 of file CompositionModel.C.

References Foam::abort(), Foam::endl(), Foam::FatalError, FatalErrorIn, forAll, phaseProperties::globalIds(), Foam::nl, phaseProperties::phase(), and WarningIn.


The documentation for this class was generated from the following files:
  • src/lagrangian/intermediate/submodels/Reacting/CompositionModel/CompositionModel/CompositionModel.H
  • src/lagrangian/intermediate/submodels/Reacting/CompositionModel/CompositionModel/CompositionModel.C
  • src/lagrangian/intermediate/submodels/Reacting/CompositionModel/CompositionModel/NewCompositionModel.C