COxidationIntrinsicRate.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) 2014-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 
27 #include "mathematicalConstants.H"
28 
29 using namespace Foam::constant;
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class CloudType>
35 (
36  const dictionary& dict,
37  CloudType& owner
38 )
39 :
40  SurfaceReactionModel<CloudType>(dict, owner, typeName),
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))),
49  CsLocalId_(-1),
50  O2GlobalId_(owner.composition().carrierId("O2")),
51  CO2GlobalId_(owner.composition().carrierId("CO2")),
52  WC_(0.0),
53  WO2_(0.0),
54  HcCO2_(0.0)
55 {
56  // Determine Cs ids
57  label idSolid = owner.composition().idSolid();
58  CsLocalId_ = owner.composition().localId(idSolid, "C");
59 
60  // Set local copies of thermo properties
61  WO2_ = owner.composition().carrier().WiValue(O2GlobalId_);
62  const scalar WCO2 = owner.composition().carrier().WiValue(CO2GlobalId_);
63  WC_ = WCO2 - WO2_;
64 
65  HcCO2_ = owner.composition().carrier().hfiValue(CO2GlobalId_);
66 
67  if (Sb_ < 0)
68  {
70  << "Stoichiometry of reaction, Sb, must be greater than zero" << nl
71  << exit(FatalError);
72  }
73 
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;
77 }
78 
79 
80 template<class CloudType>
82 (
84 )
85 :
87  Sb_(srm.Sb_),
88  C1_(srm.C1_),
89  rMean_(srm.rMean_),
90  theta_(srm.theta_),
91  Ai_(srm.Ai_),
92  Ei_(srm.Ei_),
93  Ag_(srm.Ag_),
94  tau_(srm.tau_),
95  CsLocalId_(srm.CsLocalId_),
96  O2GlobalId_(srm.O2GlobalId_),
97  CO2GlobalId_(srm.CO2GlobalId_),
98  WC_(srm.WC_),
99  WO2_(srm.WO2_),
100  HcCO2_(srm.HcCO2_)
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
105 
106 template<class CloudType>
109 {}
110 
111 
112 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
113 
114 template<class CloudType>
116 (
117  const scalar dt,
118  const label celli,
119  const scalar d,
120  const scalar T,
121  const scalar Tc,
122  const scalar pc,
123  const scalar rhoc,
124  const scalar mass,
125  const scalarField& YGas,
126  const scalarField& YLiquid,
127  const scalarField& YSolid,
128  const scalarField& YMixture,
129  const scalar N,
130  scalarField& dMassGas,
131  scalarField& dMassLiquid,
132  scalarField& dMassSolid,
133  scalarField& dMassSRCarrier
134 ) const
135 {
136  // Fraction of remaining combustible material
137  const label idSolid = CloudType::parcelType::SLD;
138  const scalar Ychar = YMixture[idSolid]*YSolid[CsLocalId_];
139 
140  // Surface combustion until combustible fraction is consumed
141  if (Ychar < small)
142  {
143  return 0.0;
144  }
145 
146  const parcelThermo& thermo = this->owner().thermo();
147  const fluidMulticomponentThermo& carrierThermo =
148  this->owner().composition().carrier();
149 
150  // Local mass fraction of O2 in the carrier phase []
151  const scalar YO2 = carrierThermo.Y(O2GlobalId_)[celli];
152 
153  // Quick exit if oxidant not present
154  if (YO2 < rootVSmall)
155  {
156  return 0.0;
157  }
158 
159  // Diffusion rate coefficient [m^2/s]
160  const scalar D0 = C1_/d*pow(0.5*(T + Tc), 0.75);
161 
162  // Apparent density of pyrolysis char [kg/m^3]
163  const scalar rhop = 6.0*mass/(constant::mathematical::pi*pow3(d));
164 
165  // Knusden diffusion coefficient [m^2/s]
166  const scalar Dkn = 97.0*rMean_*sqrt(T/WO2_);
167 
168  // Effective diffusion [m^2/s]
169  const scalar De = theta_/sqr(tau_)/(1.0/Dkn + 1/D0);
170 
171  // Cell carrier phase O2 species density [kg/m^3]
172  const scalar rhoO2 = rhoc*YO2;
173 
174  // Partial pressure O2 [Pa]
175  const scalar ppO2 = rhoO2/WO2_*RR*Tc;
176 
177  // Intrinsic reactivity [1/s]
178  const scalar ki = Ai_*exp(-Ei_/RR/T);
179 
180  // Thiele modulus []
181  const scalar phi =
182  max(0.5*d*sqrt(Sb_*rhop*Ag_*ki*ppO2/(De*rhoO2)), rootVSmall);
183 
184  // Effectiveness factor []
185  const scalar eta = max(3.0*sqr(phi)*(phi/tanh(phi) - 1.0), 0.0);
186 
187  // Chemical rate [kmol/m2/s]
188  const scalar R = eta*d/6.0*rhop*Ag_*ki;
189 
190  // Particle surface area [m^2]
191  const scalar Ap = constant::mathematical::pi*sqr(d);
192 
193  // Change in C mass [kg]
194  scalar dmC = Ap*rhoc*RR*Tc*YO2/WO2_*D0*R/(D0 + R)*dt;
195 
196  // Limit mass transfer by availability of C
197  dmC = min(mass*Ychar, dmC);
198 
199  // Molar consumption [kmol]
200  const scalar dOmega = dmC/WC_;
201 
202  // Change in O2 mass [kg]
203  const scalar dmO2 = dOmega*Sb_*WO2_;
204 
205  // Mass of newly created CO2 [kg]
206  const scalar dmCO2 = dOmega*(WC_ + Sb_*WO2_);
207 
208  // Update local particle C mass
209  dMassSolid[CsLocalId_] += dOmega*WC_;
210 
211  // Update carrier O2 and CO2 mass
212  dMassSRCarrier[O2GlobalId_] -= dmO2;
213  dMassSRCarrier[CO2GlobalId_] += dmCO2;
214 
215  const scalar hsC = thermo.solids().properties()[CsLocalId_].hs(T);
216 
217  // carrier sensible enthalpy exchange handled via change in mass
218 
219  // Heat of reaction [J]
220  return dmC*hsC - dmCO2*HcCO2_;
221 }
222 
223 
224 // ************************************************************************* //
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.
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
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:334
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
Collection of constants.
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
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
dimensionedScalar pow3(const dimensionedScalar &ds)
dimensionedScalar tanh(const dimensionedScalar &ds)
messageStream Info
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)
Definition: invIncGamma.C:102
error FatalError
static const char nl
Definition: Ostream.H:266
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dictionary dict
fluidMulticomponentThermo & thermo
Definition: createFields.H:31