33template<
class CloudType>
42 epsilon_(this->
coeffDict().getScalar(
"epsilon")),
43 gamma_(this->
coeffDict().getScalar(
"gamma")),
44 sigma_(this->
coeffDict().getScalar(
"sigma")),
47 Aeff_(this->
coeffDict().getScalar(
"Aeff")),
49 nuFuel_(this->
coeffDict().getScalar(
"nuFuel")),
50 nuOx_(this->
coeffDict().getScalar(
"nuOx")),
51 nuProd_(this->
coeffDict().getScalar(
"nuProd")),
58 label idSolid =
owner.composition().idSolid();
60 owner.composition().localId
67 owner.composition().localId
74 WO2_ =
owner.thermo().carrier().W(O2GlobalId_);
78template<
class CloudType>
86 epsilon_(srm.epsilon_),
96 O2GlobalId_(srm.O2GlobalId_),
97 FuelLocalId_(srm.FuelLocalId_),
98 ProdLocalId_(srm.ProdLocalId_),
105template<
class CloudType>
127 const scalar fComb = YSolid[FuelLocalId_];
135 const SLGThermo&
thermo = this->owner().thermo();
136 const auto&
composition = this->owner().composition();
138 const scalar WFuel =
composition.solids().properties()[FuelLocalId_].W();
139 const scalar WProd =
composition.solids().properties()[ProdLocalId_].W();
143 thermo.carrier().Y(O2GlobalId_)[celli]*rhoc/WO2_;
151 const scalar
k = A_*
exp(-Ea_/(RR*
T));
154 const scalar Deff = D12_*epsilon_/gamma_;
157 const scalar Sc =
nu/(D12_ + ROOTVSMALL);
163 const scalar r = d/2;
165 const scalar
f =
F[FuelLocalId_];
169 const scalar deltaRho0 = (nuOx_/nuFuel_)*
rhof/WFuel;
181 F[FuelLocalId_] += dfdt*dt;
184 const scalar ri = r*
cbrt(1-
f);
188 const scalar dridt = -dfdt*(
pow3(r)/3)/
sqr(ri);
194 const scalar dOmega = q02*dt;
198 composition.solids().properties()[ProdLocalId_].Hf()
199 -
composition.solids().properties()[FuelLocalId_].Hf();
202 const scalar sFuel = nuFuel_/(nuOx_);
205 const scalar sProd = nuProd_/(nuOx_);
208 dMassSRCarrier[O2GlobalId_] += dOmega*WO2_;
211 dMassSolid[FuelLocalId_] -= dOmega*WFuel*sFuel;
214 dMassSolid[ProdLocalId_] += dOmega*WProd*sProd;
218 Pout<<
"mass = " << mass <<
nl
219 <<
"fComb = " << fComb <<
nl
220 <<
"dfdt = " << dfdt <<
nl
221 <<
"F = " <<
F[FuelLocalId_] <<
nl
222 <<
"ri = " << ri <<
nl
223 <<
"dridt = " << dridt <<
nl
224 <<
"q02 = " << q02 <<
nl
225 <<
"dOmega = " << dOmega <<
nl
226 <<
"Hr = " << dOmega*WFuel*sFuel*Hc <<
endl;
230 return -dOmega*WFuel*sFuel*Hc;
234template<
class CloudType>
CGAL::indexedCell< K > Cb
const CloudType & owner() const
Return const access to the owner cloud.
Base class for heterogeneous reacting models.
HeterogeneousReactingModel(CloudType &owner)
Construct null from owner.
Heteregeneous noncatalytic reaction MUCS approach. Reference: D. Papanastassiou and G....
virtual label nReactions() const
Number of reactions in the model.
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 &YSolid, scalarField &F, const scalar N, scalar &NCpW, scalarField &dMassSolid, scalarField &dMassSRCarrier) const
Update surface reactions.
MUCSheterogeneousRate(const dictionary &dict, CloudType &owner)
Construct from dictionary.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
Density [kg/m3].
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
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
volVectorField F(fluid.F())
constexpr scalar pi(M_PI)
Namespace for handling debugging switches.
DSMCCloud< dsmcParcel > CloudType
dimensionedScalar exp(const dimensionedScalar &ds)
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
Field< scalar > scalarField
Specialisation of Field<T> for scalar.
const word GlobalIOList< Tuple2< scalar, vector > >::typeName("scalarVectorTable")
Ostream & endl(Ostream &os)
Add newline and flush stream.
dimensionedScalar sqrt(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cmplx)
Extract real component.
dimensionedScalar cbrt(const dimensionedScalar &ds)
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).
const Vector< label > N(dict.get< Vector< label > >("N"))