34template<
class CloudType>
35Foam::label Foam::COxidationMurphyShaddix<CloudType>::maxIters_ = 1000;
37template<
class CloudType>
38Foam::scalar Foam::COxidationMurphyShaddix<CloudType>::tolerance_ = 1
e-06;
43template<
class CloudType>
52 rho0_(this->
coeffDict().getScalar(
"rho0")),
58 WVol_(this->
coeffDict().getScalar(
"WVol")),
67 label idSolid =
owner.composition().idSolid();
68 CsLocalId_ =
owner.composition().localId(idSolid,
"C");
71 WO2_ =
owner.thermo().carrier().W(O2GlobalId_);
72 const scalar WCO2 =
owner.thermo().carrier().W(CO2GlobalId_);
75 HcCO2_ =
owner.thermo().carrier().Hc(CO2GlobalId_);
77 const scalar YCloc =
owner.composition().
Y0(idSolid)[CsLocalId_];
78 const scalar YSolidTot =
owner.composition().YMixture0()[idSolid];
79 Info<<
" C(s): particle mass fraction = " << YCloc*YSolidTot <<
endl;
83template<
class CloudType>
98 CsLocalId_(srm.CsLocalId_),
99 O2GlobalId_(srm.O2GlobalId_),
100 CO2GlobalId_(srm.CO2GlobalId_),
109template<
class CloudType>
134 const label idSolid = CloudType::parcelType::SLD;
135 const scalar fComb = YMixture[idSolid]*YSolid[CsLocalId_];
143 const SLGThermo&
thermo = this->owner().thermo();
146 const scalar rhoO2 = rhoc*
thermo.carrier().Y(O2GlobalId_)[celli];
158 const scalar
D = D0_*(rho0_/rhoc)*
pow(Tc/T0_, Dn_);
161 const scalar ppO2 = rhoO2/WO2_*
RR*Tc;
164 const scalar
C = pc/(
RR*Tc);
168 Pout<<
"mass = " << mass <<
nl
169 <<
"fComb = " << fComb <<
nl
170 <<
"Ap = " << Ap <<
nl
171 <<
"dt = " << dt <<
nl
180 const scalar qCsLim = mass*fComb/(WC_*Ap*dt);
184 Pout<<
"qCsLim = " << qCsLim <<
endl;
188 while ((
mag(qCs - qCsOld)/qCs > tolerance_) && (iter <= maxIters_))
191 const scalar PO2Surface = ppO2*
exp(-(qCs +
N)*d/(2*
C*
D));
192 qCs = A_*
exp(-E_/(RR*
T))*
pow(PO2Surface, n_);
193 qCs = (100.0*qCs + iter*qCsOld)/(100.0 + iter);
194 qCs =
min(qCs, qCsLim);
198 Pout<<
"iter = " << iter
199 <<
", qCsOld = " << qCsOld
207 if (iter > maxIters_)
210 <<
"iter limit reached (" << maxIters_ <<
")" <<
nl <<
endl;
214 scalar dOmega = qCs*Ap*dt;
217 dMassSRCarrier[O2GlobalId_] += -dOmega*WO2_;
218 dMassSRCarrier[CO2GlobalId_] += dOmega*(WC_ + WO2_);
221 dMassSolid[CsLocalId_] += dOmega*WC_;
223 const scalar HsC =
thermo.solids().properties()[CsLocalId_].Hs(
T);
228 return dOmega*(WC_*HsC - (WC_ + WO2_)*HcCO2_);
Limited to C(s) + O2 -> CO2.
virtual scalar calculate(const scalar dt, const scalar Re, const scalar nu, 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.
COxidationMurphyShaddix(const dictionary &dict, CloudType &owner)
Construct from dictionary.
Graphite solid properties.
const CloudType & owner() const
Return const access to the owner cloud.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Templated surface reaction model class.
SurfaceReactionModel(CloudType &owner)
Construct null from owner.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
const dictionary & dict() const
Return const access to the cloud dictionary.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
basicSpecieMixture & composition
#define WarningInFunction
Report a warning using Foam::Warning.
constexpr scalar pi(M_PI)
const scalar RR
Universal gas constant: default in [J/(kmol K)].
Namespace for handling debugging switches.
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
messageStream Info
Information stream (stdout output on master, null elsewhere).
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
scalarField Re(const UList< complex > &cmplx)
Extract real component.
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
constexpr char nl
The newline '\n' character (0x0a).
scalarList Y0(nSpecie, Zero)
const dimensionedScalar & D
psiReactionThermo & thermo
const Vector< label > N(dict.get< Vector< label > >("N"))