ReactingMultiphaseParcel< ParcelType > Class Template Reference

Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase. More...

Collaboration diagram for ReactingMultiphaseParcel< ParcelType >:
[legend]

Classes

class  constantProperties
 Class to hold reacting multiphase particle constant properties. More...
 
class  iNew
 Factory class to read-construct particles used for. More...
 

Public Types

typedef ParcelType::trackingData trackingData
 Use base tracking data. More...
 

Public Member Functions

 TypeName ("ReactingMultiphaseParcel")
 Runtime type information. More...
 
 AddToPropertyList (ParcelType, " nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)")
 String representation of properties. More...
 
 ReactingMultiphaseParcel (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
 Construct from mesh, position and topology. More...
 
 ReactingMultiphaseParcel (const polyMesh &mesh, const vector &position, const label celli)
 Construct from a position and a cell, searching for the rest of the. More...
 
 ReactingMultiphaseParcel (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti, const label typeId, const scalar nParticle0, const scalar d0, const scalar dTarget0, const vector &U0, const vector &f0, const vector &angularMomentum0, const vector &torque0, const scalarField &Y0, const scalarField &YGas0, const scalarField &YLiquid0, const scalarField &YSolid0, const constantProperties &constProps)
 Construct from components. More...
 
 ReactingMultiphaseParcel (const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
 Construct from Istream. More...
 
 ReactingMultiphaseParcel (const ReactingMultiphaseParcel &p)
 Construct as a copy. More...
 
 ReactingMultiphaseParcel (const ReactingMultiphaseParcel &p, const polyMesh &mesh)
 Construct as a copy. More...
 
virtual autoPtr< particleclone () const
 Construct and return a (basic particle) clone. More...
 
virtual autoPtr< particleclone (const polyMesh &mesh) const
 Construct and return a (basic particle) clone. More...
 
const scalarFieldYGas () const
 Return const access to mass fractions of gases. More...
 
const scalarFieldYLiquid () const
 Return const access to mass fractions of liquids. More...
 
const scalarFieldYSolid () const
 Return const access to mass fractions of solids. More...
 
label canCombust () const
 Return const access to the canCombust flag. More...
 
scalarFieldYGas ()
 Return access to mass fractions of gases. More...
 
scalarFieldYLiquid ()
 Return access to mass fractions of liquids. More...
 
scalarFieldYSolid ()
 Return access to mass fractions of solids. More...
 
label & canCombust ()
 Return access to the canCombust flag. More...
 
template<class TrackCloudType >
void setCellValues (TrackCloudType &cloud, trackingData &td)
 Set cell values. More...
 
template<class TrackCloudType >
void cellValueSourceCorrection (TrackCloudType &cloud, trackingData &td, const scalar dt)
 Correct cell values using latest transfer information. More...
 
template<class TrackCloudType >
void calc (TrackCloudType &cloud, trackingData &td, const scalar dt)
 Update parcel properties over the time interval. More...
 
void writeProperties (Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
 Write individual parcel properties to stream. More...
 
template<class TrackCloudType >
Foam::scalar CpEff (TrackCloudType &cloud, trackingData &td, const scalar p, const scalar T, const label idG, const label idL, const label idS) const
 
template<class TrackCloudType >
Foam::scalar HsEff (TrackCloudType &cloud, trackingData &td, const scalar p, const scalar T, const label idG, const label idL, const label idS) const
 
template<class TrackCloudType >
Foam::scalar LEff (TrackCloudType &cloud, trackingData &td, const scalar p, const scalar T, const label idG, const label idL, const label idS) const
 
template<class TrackCloudType >
Foam::scalar updatedDeltaVolume (TrackCloudType &cloud, const scalarField &dMassGas, const scalarField &dMassLiquid, const scalarField &dMassSolid, const label idG, const label idL, const label idS, const scalar p, const scalar T)
 

Static Public Member Functions

template<class CloudType , class CompositionType >
static void readFields (CloudType &c, const CompositionType &compModel)
 Read - composition supplied. More...
 
template<class CloudType >
static void readFields (CloudType &c)
 Read - no composition. More...
 
template<class CloudType , class CompositionType >
static void writeFields (const CloudType &c, const CompositionType &compModel)
 Write - composition supplied. More...
 
template<class CloudType >
static void writeFields (const CloudType &c)
 Read - no composition. More...
 
template<class CloudType >
static void readObjects (CloudType &c, const objectRegistry &obr)
 Read particle fields as objects from the obr registry. More...
 
template<class CloudType , class CompositionType >
static void readObjects (CloudType &c, const CompositionType &compModel, const objectRegistry &obr)
 Read particle fields as objects from the obr registry. More...
 
template<class CloudType >
static void writeObjects (const CloudType &c, objectRegistry &obr)
 Write particle fields as objects into the obr registry. More...
 
template<class CloudType , class CompositionType >
static void writeObjects (const CloudType &c, const CompositionType &compModel, objectRegistry &obr)
 Write particle fields as objects into the obr registry. More...
 

Static Public Attributes

static const std::size_t sizeofFields
 Size in bytes of the fields. More...
 
static const label GAS
 
static const label LIQ
 
static const label SLD
 

Protected Member Functions

template<class TrackCloudType >
scalar updatedDeltaVolume (TrackCloudType &cloud, const scalarField &dMassGas, const scalarField &dMassLiquid, const scalarField &dMassSolid, const label idG, const label idL, const label idS, const scalar p, const scalar T)
 Return change of volume due to mass exchange. More...
 
template<class TrackCloudType >
void calcDevolatilisation (TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar age, const scalar Ts, const scalar d, const scalar T, const scalar mass, const scalar mass0, const scalarField &YGasEff, const scalarField &YLiquidEff, const scalarField &YSolidEff, label &canCombust, scalarField &dMassDV, scalar &Sh, scalar &N, scalar &NCpW, scalarField &Cs) const
 Calculate Devolatilisation. More...
 
template<class TrackCloudType >
void calcSurfaceReactions (TrackCloudType &cloud, trackingData &td, const scalar dt, const scalar d, const scalar Re, const scalar nu, const scalar T, const scalar mass, const label canCombust, const scalar N, const scalarField &YMix, const scalarField &YGas, const scalarField &YLiquid, const scalarField &YSolid, scalarField &dMassSRGas, scalarField &dMassSRLiquid, scalarField &dMassSRSolid, scalarField &dMassSRCarrier, scalar &Sh, scalar &dhsTrans) const
 Calculate surface reactions. More...
 

Protected Attributes

scalarField YGas_
 Mass fractions of gases []. More...
 
scalarField YLiquid_
 Mass fractions of liquids []. More...
 
scalarField YSolid_
 Mass fractions of solids []. More...
 
label canCombust_
 Flag to identify if the particle can devolatilise and combust. More...
 

Friends

Ostreamoperator (Ostream &, const ReactingMultiphaseParcel< ParcelType > &)
 

Detailed Description

template<class ParcelType>
class Foam::ReactingMultiphaseParcel< ParcelType >

Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase.

Source files

Definition at line 55 of file ReactingMultiphaseParcel.H.

Member Typedef Documentation

◆ trackingData

typedef ParcelType::trackingData trackingData

Use base tracking data.

Definition at line 129 of file ReactingMultiphaseParcel.H.

Constructor & Destructor Documentation

◆ ReactingMultiphaseParcel() [1/6]

ReactingMultiphaseParcel ( const polyMesh mesh,
const barycentric coordinates,
const label  celli,
const label  tetFacei,
const label  tetPti 
)
inline

Construct from mesh, position and topology.

Other properties initialised as null

Definition at line 71 of file ReactingMultiphaseParcelI.H.

Referenced by ReactingMultiphaseParcel< ParcelType >::clone().

Here is the caller graph for this function:

◆ ReactingMultiphaseParcel() [2/6]

ReactingMultiphaseParcel ( const polyMesh mesh,
const vector position,
const label  celli 
)
inline

Construct from a position and a cell, searching for the rest of the.

required topology. Other properties are initialised as null.

Definition at line 89 of file ReactingMultiphaseParcelI.H.

◆ ReactingMultiphaseParcel() [3/6]

ReactingMultiphaseParcel ( const polyMesh mesh,
const barycentric coordinates,
const label  celli,
const label  tetFacei,
const label  tetPti,
const label  typeId,
const scalar  nParticle0,
const scalar  d0,
const scalar  dTarget0,
const vector U0,
const vector f0,
const vector angularMomentum0,
const vector torque0,
const scalarField Y0,
const scalarField YGas0,
const scalarField YLiquid0,
const scalarField YSolid0,
const constantProperties constProps 
)
inline

Construct from components.

Definition at line 105 of file ReactingMultiphaseParcelI.H.

◆ ReactingMultiphaseParcel() [4/6]

ReactingMultiphaseParcel ( const polyMesh mesh,
Istream is,
bool  readFields = true,
bool  newFormat = true 
)

Construct from Istream.

Definition at line 50 of file ReactingMultiphaseParcelIO.C.

References IOstream::check(), FUNCTION_NAME, Foam::readFields(), and DynamicList< T, SizeMin >::transfer().

Here is the call graph for this function:

◆ ReactingMultiphaseParcel() [5/6]

ReactingMultiphaseParcel ( const ReactingMultiphaseParcel< ParcelType > &  p)

Construct as a copy.

◆ ReactingMultiphaseParcel() [6/6]

ReactingMultiphaseParcel ( const ReactingMultiphaseParcel< ParcelType > &  p,
const polyMesh mesh 
)

Construct as a copy.

Member Function Documentation

◆ updatedDeltaVolume() [1/2]

scalar updatedDeltaVolume ( TrackCloudType &  cloud,
const scalarField dMassGas,
const scalarField dMassLiquid,
const scalarField dMassSolid,
const label  idG,
const label  idL,
const label  idS,
const scalar  p,
const scalar  T 
)
protected

Return change of volume due to mass exchange.

◆ calcDevolatilisation()

void calcDevolatilisation ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt,
const scalar  age,
const scalar  Ts,
const scalar  d,
const scalar  T,
const scalar  mass,
const scalar  mass0,
const scalarField YGasEff,
const scalarField YLiquidEff,
const scalarField YSolidEff,
label &  canCombust,
scalarField dMassDV,
scalar &  Sh,
scalar &  N,
scalar &  NCpW,
scalarField Cs 
) const
protected

Calculate Devolatilisation.

Definition at line 594 of file ReactingMultiphaseParcel.C.

References beta(), Foam::cbrt(), composition, Cp, Cs, forAll, Foam::max(), N(), Foam::pow(), Foam::constant::thermodynamic::RR, Foam::sqr(), Foam::sqrt(), Foam::sum(), and Foam::T().

Here is the call graph for this function:

◆ calcSurfaceReactions()

void calcSurfaceReactions ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt,
const scalar  d,
const scalar  Re,
const scalar  nu,
const scalar  T,
const scalar  mass,
const label  canCombust,
const scalar  N,
const scalarField YMix,
const scalarField YGas,
const scalarField YLiquid,
const scalarField YSolid,
scalarField dMassSRGas,
scalarField dMassSRLiquid,
scalarField dMassSRSolid,
scalarField dMassSRCarrier,
scalar &  Sh,
scalar &  dhsTrans 
) const
protected

Calculate surface reactions.

Definition at line 695 of file ReactingMultiphaseParcel.C.

References Foam::min(), N(), nu, Foam::Re(), Foam::sum(), and Foam::T().

Here is the call graph for this function:

◆ TypeName()

TypeName ( "ReactingMultiphaseParcel< ParcelType >"  )

Runtime type information.

◆ AddToPropertyList()

AddToPropertyList ( ParcelType  ,
" nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)"   
)

String representation of properties.

◆ clone() [1/2]

virtual autoPtr<particle> clone ( ) const
inlinevirtual

Construct and return a (basic particle) clone.

Definition at line 360 of file ReactingMultiphaseParcel.H.

References ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().

Here is the call graph for this function:

◆ clone() [2/2]

virtual autoPtr<particle> clone ( const polyMesh mesh) const
inlinevirtual

Construct and return a (basic particle) clone.

Definition at line 366 of file ReactingMultiphaseParcel.H.

References mesh, and ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().

Here is the call graph for this function:

◆ YGas() [1/2]

const Foam::scalarField & YGas ( ) const
inline

Return const access to mass fractions of gases.

Definition at line 191 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YGas_.

◆ YLiquid() [1/2]

const Foam::scalarField & YLiquid ( ) const
inline

Return const access to mass fractions of liquids.

Definition at line 199 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YLiquid_.

◆ YSolid() [1/2]

const Foam::scalarField & YSolid ( ) const
inline

Return const access to mass fractions of solids.

Definition at line 207 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YSolid_.

◆ canCombust() [1/2]

Foam::label canCombust ( ) const
inline

Return const access to the canCombust flag.

Definition at line 215 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::canCombust_.

◆ YGas() [2/2]

Foam::scalarField & YGas ( )
inline

Return access to mass fractions of gases.

Definition at line 222 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YGas_.

◆ YLiquid() [2/2]

Foam::scalarField & YLiquid ( )
inline

Return access to mass fractions of liquids.

Definition at line 229 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YLiquid_.

◆ YSolid() [2/2]

Foam::scalarField & YSolid ( )
inline

Return access to mass fractions of solids.

Definition at line 236 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YSolid_.

◆ canCombust() [2/2]

Foam::label & canCombust ( )
inline

Return access to the canCombust flag.

Definition at line 243 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::canCombust_.

◆ setCellValues()

void setCellValues ( TrackCloudType &  cloud,
trackingData td 
)

Set cell values.

Definition at line 189 of file ReactingMultiphaseParcel.C.

◆ cellValueSourceCorrection()

void cellValueSourceCorrection ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt 
)

Correct cell values using latest transfer information.

Definition at line 201 of file ReactingMultiphaseParcel.C.

◆ calc()

void calc ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt 
)

Update parcel properties over the time interval.

Definition at line 215 of file ReactingMultiphaseParcel.C.

References Foam::cbrt(), composition, Cs, forAll, Foam::constant::mathematical::pi(), Foam::pow4(), Foam::Re(), rho0, Su, Foam::sum(), T0, Foam::fieldTypes::volume, and Foam::Zero.

Here is the call graph for this function:

◆ readFields() [1/2]

void readFields ( CloudType c,
const CompositionType &  compModel 
)
static

Read - composition supplied.

Definition at line 91 of file ReactingMultiphaseParcelIO.C.

References Foam::constant::universal::c, forAll, Foam::max(), p, Foam::readFields(), and solidNames().

Here is the call graph for this function:

◆ readFields() [2/2]

void readFields ( CloudType c)
static

Read - no composition.

Definition at line 82 of file ReactingMultiphaseParcelIO.C.

References Foam::constant::universal::c, and Foam::readFields().

Here is the call graph for this function:

◆ writeFields() [1/2]

void writeFields ( const CloudType c,
const CompositionType &  compModel 
)
static

Write - composition supplied.

Definition at line 191 of file ReactingMultiphaseParcelIO.C.

References Foam::constant::universal::c, forAll, Foam::max(), p0, solidNames(), and Foam::writeFields().

Here is the call graph for this function:

◆ writeFields() [2/2]

void writeFields ( const CloudType c)
static

Read - no composition.

Definition at line 182 of file ReactingMultiphaseParcelIO.C.

References Foam::constant::universal::c, and Foam::writeFields().

Here is the call graph for this function:

◆ writeProperties()

void writeProperties ( Ostream os,
const wordRes filters,
const word delim,
const bool  namesOnly = false 
) const

Write individual parcel properties to stream.

Definition at line 281 of file ReactingMultiphaseParcelIO.C.

References writeProp.

◆ readObjects() [1/2]

void readObjects ( CloudType c,
const objectRegistry obr 
)
static

Read particle fields as objects from the obr registry.

  • no composition

Definition at line 306 of file ReactingMultiphaseParcelIO.C.

References Foam::constant::universal::c.

◆ readObjects() [2/2]

void readObjects ( CloudType c,
const CompositionType &  compModel,
const objectRegistry obr 
)
static

Read particle fields as objects from the obr registry.

Definition at line 330 of file ReactingMultiphaseParcelIO.C.

References Foam::constant::universal::c, forAll, Foam::max(), p0, and solidNames().

Here is the call graph for this function:

◆ writeObjects() [1/2]

void writeObjects ( const CloudType c,
objectRegistry obr 
)
static

Write particle fields as objects into the obr registry.

  • no composition

Definition at line 318 of file ReactingMultiphaseParcelIO.C.

References Foam::constant::universal::c.

◆ writeObjects() [2/2]

void writeObjects ( const CloudType c,
const CompositionType &  compModel,
objectRegistry obr 
)
static

Write particle fields as objects into the obr registry.

Definition at line 396 of file ReactingMultiphaseParcelIO.C.

References Foam::constant::universal::c, forAll, Foam::max(), p0, and solidNames().

Here is the call graph for this function:

◆ CpEff()

Foam::scalar CpEff ( TrackCloudType &  cloud,
trackingData td,
const scalar  p,
const scalar  T,
const label  idG,
const label  idL,
const label  idS 
) const

Definition at line 51 of file ReactingMultiphaseParcel.C.

References p, and Foam::T().

Here is the call graph for this function:

◆ HsEff()

Foam::scalar HsEff ( TrackCloudType &  cloud,
trackingData td,
const scalar  p,
const scalar  T,
const label  idG,
const label  idL,
const label  idS 
) const

Definition at line 71 of file ReactingMultiphaseParcel.C.

References p, and Foam::T().

Here is the call graph for this function:

◆ LEff()

Foam::scalar LEff ( TrackCloudType &  cloud,
trackingData td,
const scalar  p,
const scalar  T,
const label  idG,
const label  idL,
const label  idS 
) const

Definition at line 91 of file ReactingMultiphaseParcel.C.

References p, and Foam::T().

Here is the call graph for this function:

◆ updatedDeltaVolume() [2/2]

Foam::scalar updatedDeltaVolume ( TrackCloudType &  cloud,
const scalarField dMassGas,
const scalarField dMassLiquid,
const scalarField dMassSolid,
const label  idG,
const label  idL,
const label  idS,
const scalar  p,
const scalar  T 
)

Definition at line 145 of file ReactingMultiphaseParcel.C.

References forAll, p, Foam::sum(), Foam::T(), and Foam::Zero.

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator

Ostream& operator ( Ostream ,
const ReactingMultiphaseParcel< ParcelType > &   
)
friend

Member Data Documentation

◆ sizeofFields

const std::size_t sizeofFields
static

Size in bytes of the fields.

Definition at line 72 of file ReactingMultiphaseParcel.H.

◆ GAS

const Foam::label GAS
static

Definition at line 77 of file ReactingMultiphaseParcel.H.

◆ LIQ

const Foam::label LIQ
static

Definition at line 78 of file ReactingMultiphaseParcel.H.

◆ SLD

const Foam::label SLD
static

Definition at line 79 of file ReactingMultiphaseParcel.H.

◆ YGas_

scalarField YGas_
protected

Mass fractions of gases [].

Definition at line 192 of file ReactingMultiphaseParcel.H.

Referenced by ReactingMultiphaseParcel< ParcelType >::YGas().

◆ YLiquid_

scalarField YLiquid_
protected

Mass fractions of liquids [].

Definition at line 195 of file ReactingMultiphaseParcel.H.

Referenced by ReactingMultiphaseParcel< ParcelType >::YLiquid().

◆ YSolid_

scalarField YSolid_
protected

Mass fractions of solids [].

Definition at line 198 of file ReactingMultiphaseParcel.H.

Referenced by ReactingMultiphaseParcel< ParcelType >::YSolid().

◆ canCombust_

label canCombust_
protected

Flag to identify if the particle can devolatilise and combust.

Combustion possible only after volatile content falls below threshold value. States include: 0 = can devolatilise, cannot combust but can change 1 = can devolatilise, can combust -1 = cannot devolatilise or combust, and cannot change

Definition at line 206 of file ReactingMultiphaseParcel.H.

Referenced by ReactingMultiphaseParcel< ParcelType >::canCombust().


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