COxidationKineticDiffusionLimitedRate.C
Go to the documentation of this file.
1 /*---------------------------------------------------------------------------*\
2  ========= |
3  \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4  \\ / O peration |
5  \\ / A nd | www.openfoam.com
6  \\/ M anipulation |
7 -------------------------------------------------------------------------------
8  Copyright (C) 2011-2016 OpenFOAM Foundation
9  Copyright (C) 2019 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
30 #include "mathematicalConstants.H"
31 
32 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
33 
34 template<class CloudType>
37 (
38  const dictionary& dict,
39  CloudType& owner
40 )
41 :
42  SurfaceReactionModel<CloudType>(dict, owner, typeName),
43  Sb_(this->coeffDict().getScalar("Sb")),
44  C1_(this->coeffDict().getScalar("C1")),
45  C2_(this->coeffDict().getScalar("C2")),
46  E_(this->coeffDict().getScalar("E")),
47  CsLocalId_(-1),
48  O2GlobalId_(owner.composition().carrierId("O2")),
49  CO2GlobalId_(owner.composition().carrierId("CO2")),
50  WC_(0.0),
51  WO2_(0.0),
52  HcCO2_(0.0)
53 {
54  // Determine Cs ids
55  label idSolid = owner.composition().idSolid();
56  CsLocalId_ = owner.composition().localId(idSolid, "C");
57 
58  // Set local copies of thermo properties
59  WO2_ = owner.thermo().carrier().W(O2GlobalId_);
60  const scalar WCO2 = owner.thermo().carrier().W(CO2GlobalId_);
61  WC_ = WCO2 - WO2_;
62 
63  HcCO2_ = owner.thermo().carrier().Hc(CO2GlobalId_);
64 
65  const scalar YCloc = owner.composition().Y0(idSolid)[CsLocalId_];
66  const scalar YSolidTot = owner.composition().YMixture0()[idSolid];
67  Info<< " C(s): particle mass fraction = " << YCloc*YSolidTot << endl;
68 }
69 
70 
71 template<class CloudType>
74 (
76 )
77 :
79  Sb_(srm.Sb_),
80  C1_(srm.C1_),
81  C2_(srm.C2_),
82  E_(srm.E_),
83  CsLocalId_(srm.CsLocalId_),
84  O2GlobalId_(srm.O2GlobalId_),
85  CO2GlobalId_(srm.CO2GlobalId_),
86  WC_(srm.WC_),
87  WO2_(srm.WO2_),
88  HcCO2_(srm.HcCO2_)
89 {}
90 
91 
92 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
93 
94 template<class CloudType>
96 (
97  const scalar dt,
98  const scalar Re,
99  const scalar nu,
100  const label celli,
101  const scalar d,
102  const scalar T,
103  const scalar Tc,
104  const scalar pc,
105  const scalar rhoc,
106  const scalar mass,
107  const scalarField& YGas,
108  const scalarField& YLiquid,
109  const scalarField& YSolid,
110  const scalarField& YMixture,
111  const scalar N,
112  scalarField& dMassGas,
113  scalarField& dMassLiquid,
114  scalarField& dMassSolid,
115  scalarField& dMassSRCarrier
116 ) const
117 {
118  // Fraction of remaining combustible material
119  const label idSolid = CloudType::parcelType::SLD;
120  const scalar fComb = YMixture[idSolid]*YSolid[CsLocalId_];
121 
122  // Surface combustion active combustible fraction is consumed
123  if (fComb < SMALL)
124  {
125  return 0.0;
126  }
127 
128  const SLGThermo& thermo = this->owner().thermo();
129 
130  // Local mass fraction of O2 in the carrier phase
131  const scalar YO2 = thermo.carrier().Y(O2GlobalId_)[celli];
132 
133  // Diffusion rate coefficient
134  const scalar D0 = C1_/d*pow(0.5*(T + Tc), 0.75);
135 
136  // Kinetic rate
137  const scalar Rk = C2_*exp(-E_/(RR*Tc));
138 
139  // Particle surface area
140  const scalar Ap = constant::mathematical::pi*sqr(d);
141 
142  // Change in C mass [kg]
143  scalar dmC = Ap*rhoc*RR*Tc*YO2/WO2_*D0*Rk/(D0 + Rk)*dt;
144 
145  // Limit mass transfer by availability of C
146  dmC = min(mass*fComb, dmC);
147 
148  // Molar consumption
149  const scalar dOmega = dmC/WC_;
150 
151  // Change in O2 mass [kg]
152  const scalar dmO2 = dOmega*Sb_*WO2_;
153 
154  // Mass of newly created CO2 [kg]
155  const scalar dmCO2 = dOmega*(WC_ + Sb_*WO2_);
156 
157  // Update local particle C mass
158  dMassSolid[CsLocalId_] += dOmega*WC_;
159 
160  // Update carrier O2 and CO2 mass
161  dMassSRCarrier[O2GlobalId_] -= dmO2;
162  dMassSRCarrier[CO2GlobalId_] += dmCO2;
163 
164  const scalar HsC = thermo.solids().properties()[CsLocalId_].Hs(T);
165 
166  // carrier sensible enthalpy exchange handled via change in mass
167 
168  // Heat of reaction [J]
169  return dmC*HsC - dmCO2*HcCO2_;
170 }
171 
172 
173 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant: default in [J/(kmol K)].
Definition: thermodynamicConstants.C:46
Foam::SurfaceReactionModel
Templated surface reaction model class.
Definition: ReactingMultiphaseCloud.H:61
mathematicalConstants.H
Foam::SLGThermo
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:64
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::COxidationKineticDiffusionLimitedRate
Kinetic/diffusion limited rate surface reaction model for coal parcels. Limited to:
Definition: COxidationKineticDiffusionLimitedRate.H:54
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
Foam::COxidationKineticDiffusionLimitedRate::calculate
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.
Definition: COxidationKineticDiffusionLimitedRate.C:96
nu
volScalarField & nu
Definition: readMechanicalProperties.H:176
Foam::COxidationKineticDiffusionLimitedRate::COxidationKineticDiffusionLimitedRate
COxidationKineticDiffusionLimitedRate(const dictionary &dict, CloudType &owner)
Construct from dictionary.
Definition: COxidationKineticDiffusionLimitedRate.C:37
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::DSMCCloud
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:71
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
COxidationKineticDiffusionLimitedRate.H
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:51
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::Re
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
N
const Vector< label > N(dict.get< Vector< label >>("N"))