31 template<
class CloudType>
32 Foam::label Foam::COxidationMurphyShaddix<CloudType>::maxIters_ = 1000;
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")),
64 label idSolid = owner.composition().idSolid();
65 CsLocalId_ = owner.composition().localId(idSolid,
"C");
68 WO2_ = owner.composition().carrier().Wi(O2GlobalId_);
69 const scalar WCO2 = owner.composition().carrier().Wi(CO2GlobalId_);
72 HcCO2_ = owner.composition().carrier().Hf(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_];
149 const scalar rhoO2 = rhoc*carrier.
Y(O2GlobalId_)[celli];
161 const scalar D = D0_*(rho0_/
rhoc)*
pow(Tc/T0_, Dn_);
164 const scalar ppO2 = rhoO2/WO2_*
RR*Tc;
167 const scalar
C = pc/(
RR*Tc);
171 Pout<<
"mass = " << mass <<
nl 172 <<
"fComb = " << fComb <<
nl 173 <<
"Ap = " << Ap <<
nl 174 <<
"dt = " << dt <<
nl 183 const scalar qCsLim = mass*fComb/(WC_*Ap*dt);
187 Pout<<
"qCsLim = " << qCsLim <<
endl;
191 while ((
mag(qCs - qCsOld)/qCs > tolerance_) && (iter <= maxIters_))
194 const scalar PO2Surface = ppO2*
exp(-(qCs + N)*d/(2*C*D));
195 qCs = A_*
exp(-E_/(
RR*T))*
pow(PO2Surface, n_);
196 qCs = (100.0*qCs + iter*qCsOld)/(100.0 + iter);
197 qCs =
min(qCs, qCsLim);
201 Pout<<
"iter = " << iter
202 <<
", qCsOld = " << qCsOld
210 if (iter > maxIters_)
213 <<
"iter limit reached (" << maxIters_ <<
")" <<
nl <<
endl;
217 scalar dOmega = qCs*Ap*dt;
220 dMassSRCarrier[O2GlobalId_] += -dOmega*WO2_;
221 dMassSRCarrier[CO2GlobalId_] += dOmega*(WC_ + WO2_);
224 dMassSolid[CsLocalId_] += dOmega*WC_;
231 return dOmega*(WC_*HsC - (WC_ + WO2_)*HcCO2_);
Graphite solid properties.
fluidReactionThermo & thermo
A list of keyword definitions, which are a keyword followed by any number of values (e...
const PtrList< solidProperties > & properties() const
Return the solidProperties properties.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package, and provides:
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Specialisation of basicMixture for a mixture consisting of a number for molecular species...
Limited to C(s) + O2 -> CO2.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
dimensionedScalar exp(const dimensionedScalar &ds)
COxidationMurphyShaddix(const dictionary &dict, CloudType &owner)
Construct from dictionary.
const solidMixtureProperties & solids() const
Return reference to the global (additional) solids.
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
volScalarField rhoc(IOobject(rhocValue.name(), runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), mesh, rhocValue)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
#define WarningInFunction
Report a warning using Foam::Warning.
prefixOSstream Pout(cout, "Pout")
dimensioned< scalar > mag(const dimensioned< Type > &)
virtual ~COxidationMurphyShaddix()
Destructor.
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
Templated surface reaction model class.
Templated base class for dsmc cloud.
const dimensionedScalar e
Elementary charge.
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.