31 template<
class CloudType>
34 template<
class CloudType>
40 template<
class CloudType>
48 D0_(this->coeffDict().template lookup<scalar>(
"D0")),
49 rho0_(this->coeffDict().template lookup<scalar>(
"rho0")),
50 T0_(this->coeffDict().template lookup<scalar>(
"T0")),
51 Dn_(this->coeffDict().template lookup<scalar>(
"Dn")),
52 A_(this->coeffDict().template lookup<scalar>(
"A")),
53 E_(this->coeffDict().template lookup<scalar>(
"E")),
54 n_(this->coeffDict().template lookup<scalar>(
"n")),
55 WVol_(this->coeffDict().template lookup<scalar>(
"WVol")),
57 O2GlobalId_(owner.composition().carrierId(
"O2")),
58 CO2GlobalId_(owner.composition().carrierId(
"CO2")),
65 CsLocalId_ =
owner.composition().localId(idSolid,
"C");
68 WO2_ =
owner.composition().carrier().WiValue(O2GlobalId_);
69 const scalar WCO2 =
owner.composition().carrier().WiValue(CO2GlobalId_);
72 HcCO2_ =
owner.composition().carrier().hfiValue(CO2GlobalId_);
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>
95 CsLocalId_(srm.CsLocalId_),
96 O2GlobalId_(srm.O2GlobalId_),
97 CO2GlobalId_(srm.CO2GlobalId_),
106 template<
class CloudType>
113 template<
class CloudType>
136 const label idSolid = CloudType::parcelType::SLD;
137 const scalar fComb = YMixture[idSolid]*YSolid[CsLocalId_];
147 this->owner().composition().carrier();
150 const scalar rhoO2 = rhoc*carrierThermo.
Y(O2GlobalId_)[celli];
162 const scalar
D = D0_*(rho0_/rhoc)*
pow(Tc/T0_, Dn_);
165 const scalar ppO2 = rhoO2/WO2_*
RR*Tc;
168 const scalar
C = pc/(
RR*Tc);
172 Pout<<
"mass = " << mass <<
nl
173 <<
"fComb = " << fComb <<
nl
174 <<
"Ap = " << Ap <<
nl
175 <<
"dt = " << dt <<
nl
184 const scalar qCsLim = mass*fComb/(WC_*Ap*dt);
188 Pout<<
"qCsLim = " << qCsLim <<
endl;
192 while ((
mag(qCs - qCsOld)/qCs > tolerance_) && (iter <= maxIters_))
195 const scalar PO2Surface = ppO2*
exp(-(qCs + N)*d/(2*
C*
D));
196 qCs = A_*
exp(-E_/(
RR*
T))*
pow(PO2Surface, n_);
197 qCs = (100.0*qCs + iter*qCsOld)/(100.0 + iter);
198 qCs =
min(qCs, qCsLim);
202 Pout<<
"iter = " << iter
203 <<
", qCsOld = " << qCsOld
211 if (iter > maxIters_)
214 <<
"iter limit reached (" << maxIters_ <<
")" <<
nl <<
endl;
218 scalar dOmega = qCs*Ap*dt;
221 dMassSRCarrier[O2GlobalId_] += -dOmega*WO2_;
222 dMassSRCarrier[CO2GlobalId_] += dOmega*(WC_ + WO2_);
225 dMassSolid[CsLocalId_] += dOmega*WC_;
232 return dOmega*(WC_*hsC - (WC_ + WO2_)*HcCO2_);
static const Foam::dimensionedScalar D("D", Foam::dimTemperature, 257.14)
Limited to C(s) + O2 -> CO2.
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 ~COxidationMurphyShaddix()
Destructor.
COxidationMurphyShaddix(const dictionary &dict, CloudType &owner)
Construct from dictionary.
Graphite solid properties.
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,...
#define WarningInFunction
Report a warning using Foam::Warning.
const dimensionedScalar e
Elementary charge.
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 pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensioned< scalar > mag(const dimensioned< Type > &)
prefixOSstream Pout(cout, "Pout")
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fluidMulticomponentThermo & thermo