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


Classes | |
| class | constantProperties |
| Class to hold reacting multiphase particle constant properties. More... | |
Public Types | |
| typedef ParcelType::trackingData | trackingData |
| Use base tracking data. More... | |
Public Member Functions | |
| AddToPropertyList (ParcelType, " mass0"+" 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, const label facei) | |
| Construct from mesh, coordinates and topology. More... | |
| ReactingMultiphaseParcel (const polyMesh &mesh, const vector &position, const label celli, label &nLocateBoundaryHits) | |
| Construct from a position and a cell, searching for the rest of the. More... | |
| ReactingMultiphaseParcel (Istream &is, bool readFields=true) | |
| Construct from Istream. More... | |
| ReactingMultiphaseParcel (const ReactingMultiphaseParcel &p) | |
| Construct as a copy. More... | |
| virtual autoPtr< particle > | clone () const |
| Construct and return a clone. More... | |
| scalar | mass0 () const |
| Return const access to initial mass [kg]. More... | |
| const scalarField & | YGas () const |
| Return const access to mass fractions of gases. More... | |
| const scalarField & | YLiquid () const |
| Return const access to mass fractions of liquids. More... | |
| const scalarField & | YSolid () const |
| Return const access to mass fractions of solids. More... | |
| label | canCombust () const |
| Return const access to the canCombust flag. More... | |
| scalar & | mass0 () |
| Return access to initial mass [kg]. More... | |
| scalarField & | YGas () |
| Return access to mass fractions of gases. More... | |
| scalarField & | YLiquid () |
| Return access to mass fractions of liquids. More... | |
| scalarField & | YSolid () |
| 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... | |
| 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 | |
| static autoPtr< ReactingMultiphaseParcel > | New (Istream &is) |
| Construct from Istream and return. More... | |
| 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... | |
Protected Member Functions | |
| template<class TrackCloudType > | |
| void | calcDevolatilisation (TrackCloudType &cloud, trackingData &td, const scalar dt, 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 | |
| scalar | mass0_ |
| Initial mass [kg]. More... | |
| 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 | |
| Ostream & | operator (Ostream &, const ReactingMultiphaseParcel< ParcelType > &) |
Multiphase variant of the reacting parcel class with one/two-way coupling with the continuous phase.
Definition at line 71 of file ReactingMultiphaseParcel.H.
| typedef ParcelType::trackingData trackingData |
Use base tracking data.
Definition at line 131 of file ReactingMultiphaseParcel.H.
|
inline |
Construct from mesh, coordinates and topology.
Other properties initialised as null
Definition at line 70 of file ReactingMultiphaseParcelI.H.
Referenced by ReactingMultiphaseParcel< ParcelType >::clone(), and ReactingMultiphaseParcel< ParcelType >::New().

|
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 90 of file ReactingMultiphaseParcelI.H.
| ReactingMultiphaseParcel | ( | Istream & | is, |
| bool | readFields = true |
||
| ) |
Construct from Istream.
Definition at line 45 of file ReactingMultiphaseParcelIO.C.
References IOstream::ASCII, IOstream::check(), IOstream::format(), ReactingMultiphaseParcel< ParcelType >::mass0_, Istream::read(), ReactingMultiphaseParcel< ParcelType >::readFields(), List< T >::transfer(), ReactingMultiphaseParcel< ParcelType >::YGas_, ReactingMultiphaseParcel< ParcelType >::YLiquid_, and ReactingMultiphaseParcel< ParcelType >::YSolid_.

| ReactingMultiphaseParcel | ( | const ReactingMultiphaseParcel< ParcelType > & | p | ) |
Construct as a copy.
|
protected |
Calculate Devolatilisation.
Definition at line 501 of file ReactingMultiphaseParcel.C.
References cloud::calculate(), CompositionModel< CloudType >::carrier(), Foam::cbrt(), Cp(), multicomponentThermo::Cpi(), forAll, CompositionModel< CloudType >::idGas(), Foam::isType(), CompositionModel< CloudType >::localToCarrierId(), Foam::max(), p, Foam::pow(), Foam::constant::physicoChemical::RR, Foam::sqr(), Foam::sqrt(), Foam::sum(), Foam::T(), Foam::W(), and multicomponentThermo::WiValue().

|
protected |
Calculate surface reactions.
Definition at line 624 of file ReactingMultiphaseParcel.C.
References cloud::calculate(), Foam::isType(), Foam::min(), Foam::sum(), and Foam::T().

| AddToPropertyList | ( | ParcelType | , |
| " mass0"+" nGas(Y1..YN)"+" nLiquid(Y1..YN)"+" nSolid(Y1..YN)" | |||
| ) |
String representation of properties.
Construct and return a clone.
Definition at line 313 of file ReactingMultiphaseParcel.H.
References ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().

|
inlinestatic |
Construct from Istream and return.
Definition at line 319 of file ReactingMultiphaseParcel.H.
References ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().

|
inline |
Return const access to initial mass [kg].
Definition at line 146 of file ReactingMultiphaseParcelI.H.
|
inline |
Return const access to mass fractions of gases.
Definition at line 153 of file ReactingMultiphaseParcelI.H.
Referenced by ReactingMultiphaseParcel< ParcelType >::writeFields().

|
inline |
Return const access to mass fractions of liquids.
Definition at line 161 of file ReactingMultiphaseParcelI.H.
Referenced by ReactingMultiphaseParcel< ParcelType >::writeFields().

|
inline |
Return const access to mass fractions of solids.
Definition at line 169 of file ReactingMultiphaseParcelI.H.
Referenced by ReactingMultiphaseParcel< ParcelType >::writeFields().

|
inline |
Return const access to the canCombust flag.
Definition at line 178 of file ReactingMultiphaseParcelI.H.
|
inline |
Return access to initial mass [kg].
Definition at line 185 of file ReactingMultiphaseParcelI.H.
|
inline |
Return access to mass fractions of gases.
Definition at line 192 of file ReactingMultiphaseParcelI.H.
|
inline |
Return access to mass fractions of liquids.
Definition at line 199 of file ReactingMultiphaseParcelI.H.
|
inline |
Return access to mass fractions of solids.
Definition at line 206 of file ReactingMultiphaseParcelI.H.
|
inline |
Return access to the canCombust flag.
Definition at line 213 of file ReactingMultiphaseParcelI.H.
| void setCellValues | ( | TrackCloudType & | cloud, |
| trackingData & | td | ||
| ) |
Set cell values.
Definition at line 132 of file ReactingMultiphaseParcel.C.
| void cellValueSourceCorrection | ( | TrackCloudType & | cloud, |
| trackingData & | td, | ||
| const scalar | dt | ||
| ) |
Correct cell values using latest transfer information.
Definition at line 144 of file ReactingMultiphaseParcel.C.
| void calc | ( | TrackCloudType & | cloud, |
| trackingData & | td, | ||
| const scalar | dt | ||
| ) |
Update parcel properties over the time interval.
Definition at line 158 of file ReactingMultiphaseParcel.C.
References CompositionModel< CloudType >::carrier(), Foam::cbrt(), forAll, hs(), multicomponentThermo::hsi(), CompositionModel< CloudType >::idGas(), CompositionModel< CloudType >::idLiquid(), CompositionModel< CloudType >::idSolid(), CompositionModel< CloudType >::localToCarrierId(), Foam::constant::mathematical::pi(), Foam::pow4(), Foam::Re(), List< T >::size(), multicomponentThermo::species(), Foam::fvc::Su(), T0, and Foam::Zero.

|
static |
Read.
Definition at line 102 of file ReactingMultiphaseParcelIO.C.
References Foam::constant::universal::c, forAll, forAllIter, IOobject::MUST_READ, p, Foam::readFields(), List< T >::size(), and Foam::blendedInterfacialModel::valid().
Referenced by ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().


|
static |
Read - no composition.
Definition at line 94 of file ReactingMultiphaseParcelIO.C.
References Foam::constant::universal::c, and Foam::readFields().

|
static |
Write.
Definition at line 217 of file ReactingMultiphaseParcelIO.C.
References Foam::constant::universal::c, forAll, forAllConstIter, IOobject::NO_READ, p, regIOobject::write(), ReactingMultiphaseParcel< ParcelType >::YGas(), ReactingMultiphaseParcel< ParcelType >::YLiquid(), and ReactingMultiphaseParcel< ParcelType >::YSolid().

|
static |
Read - composition supplied.
Definition at line 209 of file ReactingMultiphaseParcelIO.C.
References Foam::constant::universal::c.
| Foam::scalar CpEff | ( | TrackCloudType & | cloud, |
| trackingData & | td, | ||
| const scalar | p, | ||
| const scalar | T, | ||
| const label | idG, | ||
| const label | idL, | ||
| const label | idS | ||
| ) | const |
| Foam::scalar hsEff | ( | TrackCloudType & | cloud, |
| trackingData & | td, | ||
| const scalar | p, | ||
| const scalar | T, | ||
| const label | idG, | ||
| const label | idL, | ||
| const label | idS | ||
| ) | const |
| Foam::scalar LEff | ( | TrackCloudType & | cloud, |
| trackingData & | td, | ||
| const scalar | p, | ||
| const scalar | T, | ||
| const label | idG, | ||
| const label | idL, | ||
| const label | idS | ||
| ) | const |
|
friend |
|
protected |
Initial mass [kg].
Definition at line 197 of file ReactingMultiphaseParcel.H.
Referenced by ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().
|
protected |
Mass fractions of gases [].
Definition at line 200 of file ReactingMultiphaseParcel.H.
Referenced by ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().
|
protected |
Mass fractions of liquids [].
Definition at line 203 of file ReactingMultiphaseParcel.H.
Referenced by ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().
|
protected |
Mass fractions of solids [].
Definition at line 206 of file ReactingMultiphaseParcel.H.
Referenced by ReactingMultiphaseParcel< ParcelType >::ReactingMultiphaseParcel().
|
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 214 of file ReactingMultiphaseParcel.H.