COxidationHurtMitchell.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration | Website: https://openfoam.org
5  \\ / A nd | Copyright (C) 2012-2023 OpenFOAM Foundation
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8 License
9  This file is part of OpenFOAM.
10 
11  OpenFOAM is free software: you can redistribute it and/or modify it
12  under the terms of the GNU General Public License as published by
13  the Free Software Foundation, either version 3 of the License, or
14  (at your option) any later version.
15 
16  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
17  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
18  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
19  for more details.
20 
21  You should have received a copy of the GNU General Public License
22  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
23 
24 \*---------------------------------------------------------------------------*/
25 
26 #include "COxidationHurtMitchell.H"
27 #include "mathematicalConstants.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class CloudType>
33 (
34  const dictionary& dict,
35  CloudType& owner
36 )
37 :
38  SurfaceReactionModel<CloudType>(dict, owner, typeName),
39  Sb_(this->coeffDict().template lookup<scalar>("Sb")),
40  CsLocalId_(-1),
41  ashLocalId_(-1),
42  O2GlobalId_(owner.composition().carrierId("O2")),
43  CO2GlobalId_(owner.composition().carrierId("CO2")),
44  WC_(0.0),
45  WO2_(0.0),
46  HcCO2_(0.0),
47  heatOfReaction_(-1.0)
48 {
49  // Determine Cs and ash ids
50  label idSolid = owner.composition().idSolid();
51  CsLocalId_ = owner.composition().localId(idSolid, "C");
52  ashLocalId_ = owner.composition().localId(idSolid, "ash", true);
53 
54  // Set local copies of thermo properties
55  WO2_ = owner.composition().carrier().WiValue(O2GlobalId_);
56  const scalar WCO2 = owner.composition().carrier().WiValue(CO2GlobalId_);
57  WC_ = WCO2 - WO2_;
58 
59  HcCO2_ = owner.composition().carrier().hfiValue(CO2GlobalId_);
60 
61  const scalar YCloc = owner.composition().Y0(idSolid)[CsLocalId_];
62  const scalar YSolidTot = owner.composition().YMixture0()[idSolid];
63  Info<< " C(s): particle mass fraction = " << YCloc*YSolidTot << endl;
64 
65  if (this->coeffDict().readIfPresent("heatOfReaction", heatOfReaction_))
66  {
67  Info<< " Using user specified heat of reaction: "
68  << heatOfReaction_ << " [J/kg]" << endl;
69  }
70 }
71 
72 
73 template<class CloudType>
75 (
77 )
78 :
80  Sb_(srm.Sb_),
81  CsLocalId_(srm.CsLocalId_),
82  ashLocalId_(srm.ashLocalId_),
83  O2GlobalId_(srm.O2GlobalId_),
84  CO2GlobalId_(srm.CO2GlobalId_),
85  WC_(srm.WC_),
86  WO2_(srm.WO2_),
87  HcCO2_(srm.HcCO2_),
88  heatOfReaction_(srm.heatOfReaction_)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
93 
94 template<class CloudType>
96 {}
97 
98 
99 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
100 
101 template<class CloudType>
103 (
104  const scalar dt,
105  const label celli,
106  const scalar d,
107  const scalar T,
108  const scalar Tc,
109  const scalar pc,
110  const scalar rhoc,
111  const scalar mass,
112  const scalarField& YGas,
113  const scalarField& YLiquid,
114  const scalarField& YSolid,
115  const scalarField& YMixture,
116  const scalar N,
117  scalarField& dMassGas,
118  scalarField& dMassLiquid,
119  scalarField& dMassSolid,
120  scalarField& dMassSRCarrier
121 ) const
122 {
123  const label idGas = CloudType::parcelType::GAS;
124  const label idSolid = CloudType::parcelType::SLD;
125  const scalar Ychar = YMixture[idSolid]*YSolid[CsLocalId_];
126 
127  // Surface combustion until combustible fraction is consumed
128  if (Ychar < small)
129  {
130  return 0.0;
131  }
132 
133  const parcelThermo& thermo = this->owner().thermo();
134  const fluidMulticomponentThermo& carrierThermo =
135  this->owner().composition().carrier();
136 
137  // Local mass fraction of O2 in the carrier phase
138  const scalar YO2 = carrierThermo.Y(O2GlobalId_)[celli];
139 
140  // No combustion if no oxygen present
141  if (YO2 < small)
142  {
143  return 0.0;
144  }
145 
146  // Conversion from [g/cm^2) to [kg/m^2]
147  const scalar convSI = 1000.0/10000.0;
148 
149  // Universal gas constant in [kcal/mol/K]
150  const scalar RRcal = 1985.877534;
151 
152  // Dry mass fraction
153  scalar Ydaf = YMixture[idGas] + YMixture[idSolid];
154  if (ashLocalId_ != -1)
155  {
156  Ydaf -= YMixture[idSolid]*YSolid[ashLocalId_];
157  }
158 
159  // Char percentage
160  const scalar charPrc = max(0, min(Ychar/(Ydaf + rootVSmall)*100.0, 100));
161 
162  // Particle surface area
163  const scalar Ap = constant::mathematical::pi*sqr(d);
164 
165  // Far field partial pressure O2 [Pa]
166  // Note: Should really use the surface partial pressure
167  const scalar ppO2 = max(0.0, rhoc*YO2/WO2_*RR*Tc);
168 
169  // Activation energy [kcal/mol]
170  const scalar E = -5.94 + 0.355*charPrc;
171 
172  // Pre-exponential factor [g/(cm^2.s.atm^0.5)]
173  const scalar lnK1750 = 2.8 - 0.0758*charPrc;
174  const scalar A = exp(lnK1750 + E/RRcal/1750.0);
175 
176  // Kinetic rate of char oxidation [g/(cm^2.s.atm^0.5)]
177  const scalar Rk = A*exp(-E/(RRcal*T));
178 
179  // Molar reaction rate per unit surface area [kmol/m^2/s]
180  const scalar qCsLim = mass*Ychar/(WC_*Ap*dt);
181  const scalar qCs = min(convSI*Rk*Foam::sqrt(ppO2/101325.0), qCsLim);
182 
183  // Calculate the number of molar units reacted [kmol]
184  const scalar dOmega = qCs*Ap*dt;
185 
186  // Add to carrier phase mass transfer
187  dMassSRCarrier[O2GlobalId_] += -dOmega*Sb_*WO2_;
188  dMassSRCarrier[CO2GlobalId_] += dOmega*(WC_ + Sb_*WO2_);
189 
190  // Add to particle mass transfer
191  dMassSolid[CsLocalId_] += dOmega*WC_;
192 
193 
194  // Return the heat of reaction [J]
195  // note: carrier sensible enthalpy exchange handled via change in mass
196  if (heatOfReaction_ < 0)
197  {
198  const scalar hsC = thermo.solids().properties()[CsLocalId_].hs(T);
199  return dOmega*(WC_*hsC - (WC_ + Sb_*WO2_)*HcCO2_);
200  }
201  else
202  {
203  return dOmega*WC_*heatOfReaction_;
204  }
205 }
206 
207 
208 // ************************************************************************* //
static const Foam::dimensionedScalar A("A", Foam::dimPressure, 611.21)
Char oxidation model given by Hurt and Mitchell:
COxidationHurtMitchell(const dictionary &dict, CloudType &owner)
Construct from dictionary.
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 ~COxidationHurtMitchell()
Destructor.
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:80
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....
Definition: dictionary.H:162
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,...
Definition: parcelThermo.H:59
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:128
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.
Definition: label.H:59
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:257
messageStream Info
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
dimensionedScalar sqrt(const dimensionedScalar &ds)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31