Loading...
Searching...
No Matches
ReactingMultiphaseParcel< ParcelType > Class Template Reference

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

#include <ReactingMultiphaseParcel.H>

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 (for parallel transfer). More...

Public Types

typedef ParcelType::trackingData trackingData
 Use base tracking data.

Public Member Functions

 TypeName ("ReactingMultiphaseParcel")
 Runtime type information.
 AddToPropertyList (ParcelType, " nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)")
 String representation of properties.
 ReactingMultiphaseParcel (const polyMesh &mesh, const barycentric &coordinates, const label celli, const label tetFacei, const label tetPti)
 Construct from mesh, position and topology.
 ReactingMultiphaseParcel (const polyMesh &mesh, const vector &position, const label celli)
 Construct from a position and a cell, searching for the rest of the.
 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.
 ReactingMultiphaseParcel (const polyMesh &mesh, Istream &is, bool readFields=true, bool newFormat=true)
 Construct from Istream.
 ReactingMultiphaseParcel (const ReactingMultiphaseParcel &p)
 Construct as a copy.
 ReactingMultiphaseParcel (const ReactingMultiphaseParcel &p, const polyMesh &mesh)
 Construct as a copy.
virtual autoPtr< particleclone () const
 Return a (basic particle) clone.
virtual autoPtr< particleclone (const polyMesh &mesh) const
 Return a (basic particle) clone.
const scalarFieldYGas () const
 Return const access to mass fractions of gases.
const scalarFieldYLiquid () const
 Return const access to mass fractions of liquids.
const scalarFieldYSolid () const
 Return const access to mass fractions of solids.
label canCombust () const
 Return const access to the canCombust flag.
scalarFieldYGas ()
 Return access to mass fractions of gases.
scalarFieldYLiquid ()
 Return access to mass fractions of liquids.
scalarFieldYSolid ()
 Return access to mass fractions of solids.
label & canCombust ()
 Return access to the canCombust flag.
template<class TrackCloudType>
void setCellValues (TrackCloudType &cloud, trackingData &td)
 Set cell values.
template<class TrackCloudType>
void cellValueSourceCorrection (TrackCloudType &cloud, trackingData &td, const scalar dt)
 Correct cell values using latest transfer information.
template<class TrackCloudType>
void calc (TrackCloudType &cloud, trackingData &td, const scalar dt)
 Update parcel properties over the time interval.
void writeProperties (Ostream &os, const wordRes &filters, const word &delim, const bool namesOnly=false) const
 Write individual parcel properties to stream.
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.
template<class CloudType>
static void readFields (CloudType &c)
 Read - no composition.
template<class CloudType, class CompositionType>
static void writeFields (const CloudType &c, const CompositionType &compModel)
 Write - composition supplied.
template<class CloudType>
static void writeFields (const CloudType &c)
 Read - no composition.
template<class CloudType>
static void readObjects (CloudType &c, const objectRegistry &obr)
 Read particle fields as objects from the obr registry.
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.
template<class CloudType>
static void writeObjects (const CloudType &c, objectRegistry &obr)
 Write particle fields as objects into the obr registry.
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.

Static Public Attributes

static const std::size_t sizeofFields
 Size in bytes of the fields.
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.
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.
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.

Protected Attributes

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

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

Member Typedef Documentation

◆ trackingData

template<class ParcelType>
typedef ParcelType::trackingData trackingData

Use base tracking data.

Definition at line 151 of file ReactingMultiphaseParcel.H.

Constructor & Destructor Documentation

◆ ReactingMultiphaseParcel() [1/6]

template<class ParcelType>
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 63 of file ReactingMultiphaseParcelI.H.

References canCombust_, coordinates(), mesh, YGas_, YLiquid_, and YSolid_.

Referenced by ReactingMultiphaseParcel< ParcelType >::iNew::operator()(), readFields(), readObjects(), writeFields(), and writeObjects().

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

◆ ReactingMultiphaseParcel() [2/6]

template<class ParcelType>
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 81 of file ReactingMultiphaseParcelI.H.

References canCombust_, mesh, YGas_, YLiquid_, and YSolid_.

◆ ReactingMultiphaseParcel() [3/6]

template<class ParcelType>
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 97 of file ReactingMultiphaseParcelI.H.

References coordinates(), mesh, and Y0().

Here is the call graph for this function:

◆ ReactingMultiphaseParcel() [4/6]

template<class ParcelType>
ReactingMultiphaseParcel ( const polyMesh & mesh,
Istream & is,
bool readFields = true,
bool newFormat = true )

Construct from Istream.

Definition at line 42 of file ReactingMultiphaseParcelIO.C.

References canCombust_, IOstream::check(), FUNCTION_NAME, mesh, readFields(), YGas_, YLiquid_, and YSolid_.

Here is the call graph for this function:

◆ ReactingMultiphaseParcel() [5/6]

template<class ParcelType>
ReactingMultiphaseParcel ( const ReactingMultiphaseParcel< ParcelType > & p)

Construct as a copy.

◆ ReactingMultiphaseParcel() [6/6]

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

Construct as a copy.

Member Function Documentation

◆ updatedDeltaVolume() [1/2]

template<class ParcelType>
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 )
protected

Return change of volume due to mass exchange.

Referenced by calc().

Here is the caller graph for this function:

◆ calcDevolatilisation()

template<class ParcelType>
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
protected

Calculate Devolatilisation.

Definition at line 586 of file ReactingMultiphaseParcel.C.

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

Referenced by calc().

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

◆ calcSurfaceReactions()

template<class ParcelType>
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
protected

Calculate surface reactions.

Definition at line 687 of file ReactingMultiphaseParcel.C.

References canCombust(), Foam::min(), N(), nu, Foam::Re(), Foam::sum(), Foam::T(), td(), YGas(), YLiquid(), and YSolid().

Referenced by calc().

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

◆ TypeName()

template<class ParcelType>
TypeName ( "ReactingMultiphaseParcel< ParcelType >" )

Runtime type information.

◆ AddToPropertyList()

template<class ParcelType>
AddToPropertyList ( ParcelType ,
" nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)"  )

String representation of properties.

◆ clone() [1/2]

template<class ParcelType>
virtual autoPtr< particle > clone ( ) const
inlinevirtual

Return a (basic particle) clone.

Definition at line 425 of file ReactingMultiphaseParcel.H.

◆ clone() [2/2]

template<class ParcelType>
virtual autoPtr< particle > clone ( const polyMesh & mesh) const
inlinevirtual

Return a (basic particle) clone.

Definition at line 433 of file ReactingMultiphaseParcel.H.

◆ YGas() [1/2]

template<class ParcelType>
const Foam::scalarField & YGas ( ) const
inline

Return const access to mass fractions of gases.

Definition at line 183 of file ReactingMultiphaseParcelI.H.

References YGas_.

Referenced by calcSurfaceReactions(), readFields(), readObjects(), writeFields(), and writeObjects().

Here is the caller graph for this function:

◆ YLiquid() [1/2]

template<class ParcelType>
const Foam::scalarField & YLiquid ( ) const
inline

Return const access to mass fractions of liquids.

Definition at line 191 of file ReactingMultiphaseParcelI.H.

References YLiquid_.

Referenced by calcSurfaceReactions(), readFields(), readObjects(), writeFields(), and writeObjects().

Here is the caller graph for this function:

◆ YSolid() [1/2]

template<class ParcelType>
const Foam::scalarField & YSolid ( ) const
inline

Return const access to mass fractions of solids.

Definition at line 199 of file ReactingMultiphaseParcelI.H.

References YSolid_.

Referenced by calcSurfaceReactions(), readFields(), readObjects(), writeFields(), and writeObjects().

Here is the caller graph for this function:

◆ canCombust() [1/2]

template<class ParcelType>
Foam::label canCombust ( ) const
inline

Return const access to the canCombust flag.

Definition at line 208 of file ReactingMultiphaseParcelI.H.

References canCombust_.

Referenced by calcDevolatilisation(), and calcSurfaceReactions().

Here is the caller graph for this function:

◆ YGas() [2/2]

template<class ParcelType>
Foam::scalarField & YGas ( )
inline

Return access to mass fractions of gases.

Definition at line 215 of file ReactingMultiphaseParcelI.H.

References YGas_.

◆ YLiquid() [2/2]

template<class ParcelType>
Foam::scalarField & YLiquid ( )
inline

Return access to mass fractions of liquids.

Definition at line 222 of file ReactingMultiphaseParcelI.H.

References YLiquid_.

◆ YSolid() [2/2]

template<class ParcelType>
Foam::scalarField & YSolid ( )
inline

Return access to mass fractions of solids.

Definition at line 229 of file ReactingMultiphaseParcelI.H.

References YSolid_.

◆ canCombust() [2/2]

template<class ParcelType>
Foam::label & canCombust ( )
inline

Return access to the canCombust flag.

Definition at line 236 of file ReactingMultiphaseParcelI.H.

References canCombust_.

◆ setCellValues()

template<class ParcelType>
template<class TrackCloudType>
void setCellValues ( TrackCloudType & cloud,
trackingData & td )

Set cell values.

Definition at line 181 of file ReactingMultiphaseParcel.C.

References td().

Here is the call graph for this function:

◆ cellValueSourceCorrection()

template<class ParcelType>
template<class TrackCloudType>
void cellValueSourceCorrection ( TrackCloudType & cloud,
trackingData & td,
const scalar dt )

Correct cell values using latest transfer information.

Definition at line 193 of file ReactingMultiphaseParcel.C.

References td().

Here is the call graph for this function:

◆ calc()

template<class ParcelType>
template<class TrackCloudType>
void calc ( TrackCloudType & cloud,
trackingData & td,
const scalar dt )

Update parcel properties over the time interval.

Definition at line 207 of file ReactingMultiphaseParcel.C.

References calcDevolatilisation(), calcSurfaceReactions(), canCombust_, Foam::cbrt(), composition, Cs, forAll, GAS, LIQ, Foam::constant::mathematical::pi(), Foam::pow4(), Foam::Re(), rho0, SLD, Foam::Su(), Foam::sum(), T0, td(), updatedDeltaVolume(), YGas_, YLiquid_, YSolid_, and Foam::Zero.

Here is the call graph for this function:

◆ readFields() [1/2]

template<class ParcelType>
template<class CloudType, class CompositionType>
void readFields ( CloudType & c,
const CompositionType & compModel )
static

Read - composition supplied.

Definition at line 83 of file ReactingMultiphaseParcelIO.C.

References forAll, GAS, LIQ, Foam::max(), IOobjectOption::MUST_READ, p, ReactingMultiphaseParcel(), UList< T >::size(), SLD, solidNames(), YGas(), YLiquid(), and YSolid().

Referenced by ReactingMultiphaseParcel().

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

◆ readFields() [2/2]

template<class ParcelType>
template<class CloudType>
void readFields ( CloudType & c)
static

Read - no composition.

Definition at line 75 of file ReactingMultiphaseParcelIO.C.

◆ writeFields() [1/2]

template<class ParcelType>
template<class CloudType, class CompositionType>
void writeFields ( const CloudType & c,
const CompositionType & compModel )
static

Write - composition supplied.

Definition at line 183 of file ReactingMultiphaseParcelIO.C.

References forAll, GAS, LIQ, Foam::max(), IOobjectOption::NO_READ, p0, ReactingMultiphaseParcel(), SLD, solidNames(), YGas(), YLiquid(), and YSolid().

Here is the call graph for this function:

◆ writeFields() [2/2]

template<class ParcelType>
template<class CloudType>
void writeFields ( const CloudType & c)
static

Read - no composition.

Definition at line 175 of file ReactingMultiphaseParcelIO.C.

◆ writeProperties()

template<class ParcelType>
void writeProperties ( Ostream & os,
const wordRes & filters,
const word & delim,
const bool namesOnly = false ) const

Write individual parcel properties to stream.

Definition at line 274 of file ReactingMultiphaseParcelIO.C.

References canCombust_, os(), writeProp, YGas_, YLiquid_, and YSolid_.

Here is the call graph for this function:

◆ readObjects() [1/2]

template<class ParcelType>
template<class CloudType>
void readObjects ( CloudType & c,
const objectRegistry & obr )
static

Read particle fields as objects from the obr registry.

  • no composition

Definition at line 299 of file ReactingMultiphaseParcelIO.C.

◆ readObjects() [2/2]

template<class ParcelType>
template<class CloudType, class CompositionType>
void readObjects ( CloudType & c,
const CompositionType & compModel,
const objectRegistry & obr )
static

Read particle fields as objects from the obr registry.

Definition at line 323 of file ReactingMultiphaseParcelIO.C.

References forAll, GAS, LIQ, cloud::lookupIOField(), Foam::max(), p0, ReactingMultiphaseParcel(), SLD, solidNames(), YGas(), YLiquid(), and YSolid().

Here is the call graph for this function:

◆ writeObjects() [1/2]

template<class ParcelType>
template<class CloudType>
void writeObjects ( const CloudType & c,
objectRegistry & obr )
static

Write particle fields as objects into the obr registry.

  • no composition

Definition at line 311 of file ReactingMultiphaseParcelIO.C.

◆ writeObjects() [2/2]

template<class ParcelType>
template<class CloudType, class CompositionType>
void writeObjects ( const CloudType & c,
const CompositionType & compModel,
objectRegistry & obr )
static

Write particle fields as objects into the obr registry.

Definition at line 390 of file ReactingMultiphaseParcelIO.C.

References cloud::createIOField(), forAll, GAS, LIQ, Foam::max(), p0, ReactingMultiphaseParcel(), SLD, solidNames(), YGas(), YLiquid(), and YSolid().

Here is the call graph for this function:

◆ CpEff()

template<class ParcelType>
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

Definition at line 43 of file ReactingMultiphaseParcel.C.

◆ HsEff()

template<class ParcelType>
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

Definition at line 63 of file ReactingMultiphaseParcel.C.

◆ LEff()

template<class ParcelType>
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

Definition at line 83 of file ReactingMultiphaseParcel.C.

◆ updatedDeltaVolume() [2/2]

template<class ParcelType>
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 )

Definition at line 137 of file ReactingMultiphaseParcel.C.

◆ operator<<

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

Member Data Documentation

◆ sizeofFields

template<class ParcelType>
const std::size_t sizeofFields
static

Size in bytes of the fields.

Definition at line 70 of file ReactingMultiphaseParcel.H.

◆ GAS

template<class ParcelType>
const Foam::label GAS
static

◆ LIQ

template<class ParcelType>
const Foam::label LIQ
static

Definition at line 76 of file ReactingMultiphaseParcel.H.

Referenced by calc(), readFields(), readObjects(), writeFields(), and writeObjects().

◆ SLD

template<class ParcelType>
const Foam::label SLD
static

Definition at line 77 of file ReactingMultiphaseParcel.H.

Referenced by calc(), readFields(), readObjects(), writeFields(), and writeObjects().

◆ YGas_

template<class ParcelType>
scalarField YGas_
protected

◆ YLiquid_

template<class ParcelType>
scalarField YLiquid_
protected

◆ YSolid_

template<class ParcelType>
scalarField YSolid_
protected

◆ canCombust_

template<class ParcelType>
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 245 of file ReactingMultiphaseParcel.H.

Referenced by calc(), canCombust(), canCombust(), ReactingMultiphaseParcel(), ReactingMultiphaseParcel(), ReactingMultiphaseParcel(), and writeProperties().


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