33 template<
class CloudType>
41 Sb_(this->coeffDict().template lookup<scalar>(
"Sb")),
42 C1_(this->coeffDict().template lookup<scalar>(
"C1")),
43 rMean_(this->coeffDict().template lookup<scalar>(
"rMean")),
44 theta_(this->coeffDict().template lookup<scalar>(
"theta")),
45 Ai_(this->coeffDict().template lookup<scalar>(
"Ai")),
46 Ei_(this->coeffDict().template lookup<scalar>(
"Ei")),
47 Ag_(this->coeffDict().template lookup<scalar>(
"Ag")),
48 tau_(this->coeffDict().lookupOrDefault(
"tau",
sqrt(2.0))),
50 O2GlobalId_(owner.composition().carrierId(
"O2")),
51 CO2GlobalId_(owner.composition().carrierId(
"CO2")),
58 CsLocalId_ =
owner.composition().localId(idSolid,
"C");
61 WO2_ =
owner.composition().carrier().WiValue(O2GlobalId_);
62 const scalar WCO2 =
owner.composition().carrier().WiValue(CO2GlobalId_);
65 HcCO2_ =
owner.composition().carrier().hfiValue(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>
95 CsLocalId_(srm.CsLocalId_),
96 O2GlobalId_(srm.O2GlobalId_),
97 CO2GlobalId_(srm.CO2GlobalId_),
106 template<
class CloudType>
114 template<
class CloudType>
137 const label idSolid = CloudType::parcelType::SLD;
138 const scalar Ychar = YMixture[idSolid]*YSolid[CsLocalId_];
148 this->owner().composition().carrier();
151 const scalar YO2 = carrierThermo.
Y(O2GlobalId_)[celli];
154 if (YO2 < rootVSmall)
160 const scalar D0 = C1_/d*
pow(0.5*(
T + Tc), 0.75);
166 const scalar Dkn = 97.0*rMean_*
sqrt(
T/WO2_);
169 const scalar De = theta_/
sqr(tau_)/(1.0/Dkn + 1/D0);
172 const scalar rhoO2 = rhoc*YO2;
175 const scalar ppO2 = rhoO2/WO2_*
RR*Tc;
178 const scalar ki = Ai_*
exp(-Ei_/
RR/
T);
182 max(0.5*d*
sqrt(Sb_*rhop*Ag_*ki*ppO2/(De*rhoO2)), rootVSmall);
185 const scalar eta =
max(3.0*
sqr(phi)*(phi/
tanh(phi) - 1.0), 0.0);
188 const scalar
R = eta*d/6.0*rhop*Ag_*ki;
194 scalar dmC = Ap*rhoc*
RR*Tc*YO2/WO2_*D0*
R/(D0 +
R)*dt;
197 dmC =
min(mass*Ychar, dmC);
200 const scalar dOmega = dmC/WC_;
203 const scalar dmO2 = dOmega*Sb_*WO2_;
206 const scalar dmCO2 = dOmega*(WC_ + Sb_*WO2_);
209 dMassSolid[CsLocalId_] += dOmega*WC_;
212 dMassSRCarrier[O2GlobalId_] -= dmO2;
213 dMassSRCarrier[CO2GlobalId_] += dmCO2;
220 return dmC*hsC - dmCO2*HcCO2_;
Intrinsic char surface reaction mndel.
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.
COxidationIntrinsicRate(const dictionary &dict, CloudType &owner)
Construct from dictionary.
virtual ~COxidationIntrinsicRate()
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,...
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
errorManipArg< error, int > exit(error &err, const int errNo=1)
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.
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
fluidMulticomponentThermo & thermo