Go to the documentation of this file.
33 template<
class CloudType>
48 tau_(this->coeffDict().lookupOrDefault(
"tau",
sqrt(2.0))),
50 O2GlobalId_(owner.composition().carrierId(
"O2")),
51 CO2GlobalId_(owner.composition().carrierId(
"CO2")),
57 label idSolid = owner.composition().idSolid();
58 CsLocalId_ = owner.composition().localId(idSolid,
"C");
61 WO2_ = owner.thermo().carrier().W(O2GlobalId_);
62 const scalar WCO2 = owner.thermo().carrier().W(CO2GlobalId_);
65 HcCO2_ = owner.thermo().carrier().Hc(CO2GlobalId_);
70 <<
"Stoichiometry of reaction, Sb, must be greater than zero" <<
nl
74 const scalar YCloc = owner.composition().Y0(idSolid)[CsLocalId_];
75 const scalar YSolidTot = owner.composition().YMixture0()[idSolid];
76 Info<<
" C(s): particle mass fraction = " << YCloc*YSolidTot <<
endl;
80 template<
class CloudType>
106 template<
class CloudType>
114 template<
class CloudType>
137 const label idSolid = CloudType::parcelType::SLD;
138 const scalar Ychar = YMixture[idSolid]*YSolid[CsLocalId_];
149 const scalar YO2 =
thermo.carrier().Y(O2GlobalId_)[cellI];
152 if (YO2 < ROOTVSMALL)
158 const scalar D0 = C1_/d*
pow(0.5*(
T + Tc), 0.75);
164 const scalar Dkn = 97.0*rMean_*
sqrt(
T/WO2_);
167 const scalar De = theta_/
sqr(tau_)/(1.0/Dkn + 1/D0);
170 const scalar rhoO2 = rhoc*YO2;
173 const scalar ppO2 = rhoO2/WO2_*
RR*Tc;
176 const scalar ki = Ai_*
exp(-Ei_/
RR/
T);
180 max(0.5*d*
sqrt(Sb_*rhop*Ag_*ki*ppO2/(De*rhoO2)), ROOTVSMALL);
186 const scalar
R = eta*d/6.0*rhop*Ag_*ki;
192 scalar dmC = Ap*rhoc*
RR*Tc*YO2/WO2_*D0*
R/(D0 +
R)*dt;
195 dmC =
min(mass*Ychar, dmC);
198 const scalar dOmega = dmC/WC_;
201 const scalar dmO2 = dOmega*Sb_*WO2_;
204 const scalar dmCO2 = dOmega*(WC_ + Sb_*WO2_);
207 dMassSolid[CsLocalId_] += dOmega*WC_;
210 dMassSRCarrier[O2GlobalId_] -= dmO2;
211 dMassSRCarrier[CO2GlobalId_] += dmCO2;
213 const scalar HsC =
thermo.solids().properties()[CsLocalId_].Hs(
T);
218 return dmC*HsC - dmCO2*HcCO2_;
const scalar RR
Universal gas constant (default in [J/(kmol K)])
Intrinsic char surface reaction mndel.
const scalar tau_
Pore tortuosity []; default to sqrt(2)
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
const scalar C1_
Mass diffusion limited rate constant.
virtual ~COxidationIntrinsicRate()
Destructor.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar exp(const dimensionedScalar &ds)
const scalar rMean_
Mean pore radius [m].
#define R(A, B, C, D, E, F, K, M)
intWM_LABEL_SIZE_t label
A label is an int32_t or int64_t as specified by the pre-processor macro WM_LABEL_SIZE.
Pre-declare SubField and related Field type.
const scalar Sb_
Stoichiometry of reaction [].
dimensionedScalar tanh(const dimensionedScalar &ds)
dimensionedScalar pow3(const dimensionedScalar &ds)
const scalar theta_
Char porosity [] = 1 - rho_apparent/rho_true.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
COxidationIntrinsicRate(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.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Templated base class for dsmc cloud.
scalar HcCO2_
Formation enthalpy for CO2 [J/kg].
A list of keyword definitions, which are a keyword followed by any number of values (e....
label CsLocalId_
Cs positions in global/local lists.
errorManipArg< error, int > exit(error &err, const int errNo=1)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
bool readScalar(const char *buf, doubleScalar &s)
Read whole of buf as a scalar. Return true if succesful.
label O2GlobalId_
O2 position in global list.
dimensionedScalar sqrt(const dimensionedScalar &ds)
const scalar Ag_
Char specific internal area [m2/kg].
label CO2GlobalId_
CO2 positions in global list.
const scalar Ai_
Pre-exponential factor.
scalar WO2_
Molecular weight of O2 [kg/kmol].
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
scalar WC_
Molecular weight of C [kg/kmol].
stressControl lookup("compactNormalStress") >> compactNormalStress
const scalar Ei_
Activation energy.