31 template<
class CloudType>
39 Sb_(this->coeffDict().template lookup<scalar>(
"Sb")),
42 O2GlobalId_(owner.composition().carrierId(
"O2")),
43 CO2GlobalId_(owner.composition().carrierId(
"CO2")),
51 CsLocalId_ =
owner.composition().localId(idSolid,
"C");
52 ashLocalId_ =
owner.composition().localId(idSolid,
"ash",
true);
55 WO2_ =
owner.composition().carrier().WiValue(O2GlobalId_);
56 const scalar WCO2 =
owner.composition().carrier().WiValue(CO2GlobalId_);
59 HcCO2_ =
owner.composition().carrier().hfiValue(CO2GlobalId_);
61 const scalar YCloc =
owner.composition().Y0(idSolid)[CsLocalId_];
62 const scalar YSolidTot =
owner.composition().YMixture0()[idSolid];
63 Info<<
" C(s): particle mass fraction = " << YCloc*YSolidTot <<
endl;
65 if (this->
coeffDict().readIfPresent(
"heatOfReaction", heatOfReaction_))
67 Info<<
" Using user specified heat of reaction: "
68 << heatOfReaction_ <<
" [J/kg]" <<
endl;
73 template<
class CloudType>
81 CsLocalId_(srm.CsLocalId_),
82 ashLocalId_(srm.ashLocalId_),
83 O2GlobalId_(srm.O2GlobalId_),
84 CO2GlobalId_(srm.CO2GlobalId_),
88 heatOfReaction_(srm.heatOfReaction_)
94 template<
class CloudType>
101 template<
class CloudType>
123 const label idGas = CloudType::parcelType::GAS;
124 const label idSolid = CloudType::parcelType::SLD;
125 const scalar Ychar = YMixture[idSolid]*YSolid[CsLocalId_];
135 this->owner().composition().carrier();
138 const scalar YO2 = carrierThermo.
Y(O2GlobalId_)[celli];
147 const scalar convSI = 1000.0/10000.0;
150 const scalar RRcal = 1985.877534;
153 scalar Ydaf = YMixture[idGas] + YMixture[idSolid];
154 if (ashLocalId_ != -1)
156 Ydaf -= YMixture[idSolid]*YSolid[ashLocalId_];
160 const scalar charPrc =
max(0,
min(Ychar/(Ydaf + rootVSmall)*100.0, 100));
167 const scalar ppO2 =
max(0.0, rhoc*YO2/WO2_*
RR*Tc);
170 const scalar E = -5.94 + 0.355*charPrc;
173 const scalar lnK1750 = 2.8 - 0.0758*charPrc;
174 const scalar
A =
exp(lnK1750 + E/RRcal/1750.0);
177 const scalar Rk =
A*
exp(-E/(RRcal*
T));
180 const scalar qCsLim = mass*Ychar/(WC_*Ap*dt);
181 const scalar qCs =
min(convSI*Rk*
Foam::sqrt(ppO2/101325.0), qCsLim);
184 const scalar dOmega = qCs*Ap*dt;
187 dMassSRCarrier[O2GlobalId_] += -dOmega*Sb_*WO2_;
188 dMassSRCarrier[CO2GlobalId_] += dOmega*(WC_ + Sb_*WO2_);
191 dMassSolid[CsLocalId_] += dOmega*WC_;
196 if (heatOfReaction_ < 0)
199 return dOmega*(WC_*hsC - (WC_ + Sb_*WO2_)*HcCO2_);
203 return dOmega*WC_*heatOfReaction_;
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
Char oxidation model given by Hurt and Mitchell:
COxidationHurtMitchell(const dictionary &dict, CloudType &owner)
Construct from dictionary.
virtual scalar calculate(const scalar dt, const label celli, const scalar d, const scalar T, const scalar Tc, const scalar pc, const scalar rhoc, const scalar mass, const scalarField &YGas, const scalarField &YLiquid, const scalarField &YSolid, const scalarField &YMixture, const scalar N, scalarField &dMassGas, scalarField &dMassLiquid, scalarField &dMassSolid, scalarField &dMassSRCarrier) const
Update surface reactions.
virtual ~COxidationHurtMitchell()
Destructor.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Templated surface reaction model class.
virtual const IOdictionary & properties() const =0
Properties dictionary.
A list of keyword definitions, which are a keyword followed by any number of values (e....
Base-class for multi-component fluid thermodynamic properties.
virtual PtrList< volScalarField > & Y()=0
Access the mass-fraction fields.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
dimensionedScalar exp(const dimensionedScalar &ds)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fluidMulticomponentThermo & thermo