ReactingMultiphaseParcel< ParcelType > Class Template Reference

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

Inheritance diagram for ReactingMultiphaseParcel< ParcelType >:
Collaboration diagram for ReactingMultiphaseParcel< ParcelType >:

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

Static Public Member Functions

template<class CloudType , class CompositionType >
static void readFields (CloudType &c, const CompositionType &compModel)
 Read. 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. More...
 
template<class CloudType >
static void writeFields (const CloudType &c)
 Read - composition supplied. More...
 

Static Public Attributes

static const label GAS
 
static const label LIQ
 
static const label SLD
 

Protected Member Functions

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 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 50 of file ReactingMultiphaseParcel.H.

Member Typedef Documentation

◆ trackingData

typedef ParcelType::trackingData trackingData

Use base tracking data.

Definition at line 130 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

◆ 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 87 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().

Here is the call graph for this function:

◆ 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 103 of file ReactingMultiphaseParcelI.H.

◆ ReactingMultiphaseParcel() [4/6]

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

Construct from Istream.

Definition at line 46 of file ReactingMultiphaseParcelIO.C.

References IOstream::check(), and DynamicList< T, SizeInc, SizeMult, SizeDiv >::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

◆ 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

◆ calcSurfaceReactions()

void calcSurfaceReactions ( TrackCloudType &  cloud,
trackingData td,
const scalar  dt,
const scalar  d,
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

◆ 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 342 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 348 of file ReactingMultiphaseParcel.H.

References 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 189 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YGas_, and ReactingMultiphaseParcel< ParcelType >::YLiquid().

Referenced by ReactingMultiphaseParcel< ParcelType >::constantProperties::hRetentionCoeff(), ReactingMultiphaseParcel< ParcelType >::iNew::operator()(), and ReactingMultiphaseParcel< ParcelType >::writeFields().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ YLiquid() [1/2]

const Foam::scalarField & YLiquid ( ) const
inline

Return const access to mass fractions of liquids.

Definition at line 197 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YLiquid_, and ReactingMultiphaseParcel< ParcelType >::YSolid().

Referenced by ReactingMultiphaseParcel< ParcelType >::iNew::operator()(), ReactingMultiphaseParcel< ParcelType >::writeFields(), and ReactingMultiphaseParcel< ParcelType >::YGas().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ YSolid() [1/2]

const Foam::scalarField & YSolid ( ) const
inline

Return const access to mass fractions of solids.

Definition at line 205 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YSolid_.

Referenced by ReactingMultiphaseParcel< ParcelType >::iNew::operator()(), ReactingMultiphaseParcel< ParcelType >::writeFields(), and ReactingMultiphaseParcel< ParcelType >::YLiquid().

Here is the caller graph for this function:

◆ canCombust() [1/2]

Foam::label canCombust ( ) const
inline

Return const access to the canCombust flag.

Definition at line 213 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::canCombust_.

Referenced by ReactingMultiphaseParcel< ParcelType >::iNew::operator()().

Here is the caller graph for this function:

◆ YGas() [2/2]

Foam::scalarField & YGas ( )
inline

Return access to mass fractions of gases.

Definition at line 220 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 227 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 234 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::YSolid_.

◆ canCombust() [2/2]

Foam::label & canCombust ( )
inline

Return access to the canCombust flag.

Definition at line 241 of file ReactingMultiphaseParcelI.H.

References ReactingMultiphaseParcel< ParcelType >::canCombust_.

◆ setCellValues()

void setCellValues ( TrackCloudType &  cloud,
trackingData td 
)

Set cell values.

Definition at line 138 of file ReactingMultiphaseParcel.C.

References ReactingMultiphaseParcel< ParcelType >::cellValueSourceCorrection().

Referenced by ReactingMultiphaseParcel< ParcelType >::LEff(), and ReactingMultiphaseParcel< ParcelType >::iNew::operator()().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ cellValueSourceCorrection()

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

Correct cell values using latest transfer information.

Definition at line 150 of file ReactingMultiphaseParcel.C.

References ReactingMultiphaseParcel< ParcelType >::calc().

Referenced by ReactingMultiphaseParcel< ParcelType >::iNew::operator()(), and ReactingMultiphaseParcel< ParcelType >::setCellValues().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ calc()

◆ readFields() [1/2]

◆ readFields() [2/2]

void readFields ( CloudType c)
static

Read - no composition.

Definition at line 92 of file ReactingMultiphaseParcelIO.C.

References Foam::readFields(), and ReactingMultiphaseParcel< ParcelType >::readFields().

Here is the call graph for this function:

◆ writeFields() [1/2]

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

◆ writeFields() [2/2]

void writeFields ( const CloudType c)
static

Read - composition supplied.

Definition at line 208 of file ReactingMultiphaseParcelIO.C.

References ReactingMultiphaseParcel< ParcelType >::writeFields().

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 48 of file ReactingMultiphaseParcel.C.

◆ 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 68 of file ReactingMultiphaseParcel.C.

◆ 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 88 of file ReactingMultiphaseParcel.C.

References Foam::max(), and ReactingMultiphaseParcel< ParcelType >::setCellValues().

Here is the call graph for this function:

Friends And Related Function Documentation

◆ operator

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

Member Data Documentation

◆ GAS

const Foam::label GAS
static

Definition at line 78 of file ReactingMultiphaseParcel.H.

◆ LIQ

const Foam::label LIQ
static

Definition at line 79 of file ReactingMultiphaseParcel.H.

◆ SLD

const Foam::label SLD
static

Definition at line 80 of file ReactingMultiphaseParcel.H.

◆ YGas_

◆ YLiquid_

◆ 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 207 of file ReactingMultiphaseParcel.H.

Referenced by ReactingMultiphaseParcel< ParcelType >::calcSurfaceReactions(), and ReactingMultiphaseParcel< ParcelType >::canCombust().


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