MUCSheterogeneousRate.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) 2018-2019 OpenCFD Ltd.
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CloudType>
35(
36 const dictionary& dict,
37 CloudType& owner
38)
39:
40 HeterogeneousReactingModel<CloudType>(dict, owner, typeName),
41 D12_(this->coeffDict().getScalar("D12")),
42 epsilon_(this->coeffDict().getScalar("epsilon")),
43 gamma_(this->coeffDict().getScalar("gamma")),
44 sigma_(this->coeffDict().getScalar("sigma")),
45 E_(this->coeffDict().getScalar("E")),
46 A_(this->coeffDict().getScalar("A")),
47 Aeff_(this->coeffDict().getScalar("Aeff")),
48 Ea_(this->coeffDict().getScalar("Ea")),
49 nuFuel_(this->coeffDict().getScalar("nuFuel")),
50 nuOx_(this->coeffDict().getScalar("nuOx")),
51 nuProd_(this->coeffDict().getScalar("nuProd")),
52 O2GlobalId_(owner.composition().carrierId("O2")),
53 FuelLocalId_(-1),
54 ProdLocalId_(-1),
55 WO2_(0.0)
56{
57 // Determine Cs ids
58 label idSolid = owner.composition().idSolid();
59 FuelLocalId_ =
60 owner.composition().localId
61 (
62 idSolid,
63 this->coeffDict().getWord("fuel")
64 );
65
66 ProdLocalId_ =
67 owner.composition().localId
68 (
69 idSolid,
70 this->coeffDict().getWord("product")
71 );
72
73 // Set local copies of thermo properties
74 WO2_ = owner.thermo().carrier().W(O2GlobalId_);
75}
76
77
78template<class CloudType>
80(
82)
83:
85 D12_(srm.D12_),
86 epsilon_(srm.epsilon_),
87 gamma_(srm.gamma_),
88 sigma_(srm.sigma_),
89 E_(srm.E_),
90 A_(srm.A_),
91 Aeff_(srm.Aeff_),
92 Ea_(srm.Ea_),
93 nuFuel_(srm.nuFuel_),
94 nuOx_(srm.nuOx_),
95 nuProd_(srm.nuProd_),
96 O2GlobalId_(srm.O2GlobalId_),
97 FuelLocalId_(srm.FuelLocalId_),
98 ProdLocalId_(srm.ProdLocalId_),
99 WO2_(srm.WO2_)
100{}
101
102
103// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
104
105template<class CloudType>
107(
108 const scalar dt,
109 const scalar Re,
110 const scalar nu,
111 const label celli,
112 const scalar d,
113 const scalar T,
114 const scalar Tc,
115 const scalar pc,
116 const scalar rhoc,
117 const scalar mass,
118 const scalarField& YSolid,
119 scalarField& F,
120 const scalar N,
121 scalar& NCpW,
122 scalarField& dMassSolid,
123 scalarField& dMassSRCarrier
124) const
125{
126 // Fraction of remaining combustible material
127 const scalar fComb = YSolid[FuelLocalId_];
128
129 // Surface combustion until combustible fraction is consumed
130 if (fComb < SMALL)
131 {
132 return 0.0;
133 }
134
135 const SLGThermo& thermo = this->owner().thermo();
136 const auto& composition = this->owner().composition();
137
138 const scalar WFuel = composition.solids().properties()[FuelLocalId_].W();
139 const scalar WProd = composition.solids().properties()[ProdLocalId_].W();
140
141 // O2 concentration [Kmol/m3]
142 const scalar Cb =
143 thermo.carrier().Y(O2GlobalId_)[celli]*rhoc/WO2_;
144
145 if (Cb < SMALL)
146 {
147 return 0.0;
148 }
149
150 // Reaction constant
151 const scalar k = A_*exp(-Ea_/(RR*T));
152
153 // Effective diffussivity
154 const scalar Deff = D12_*epsilon_/gamma_;
155
156 // Schmidt number
157 const scalar Sc = nu/(D12_ + ROOTVSMALL);
158
159 // Mass transfer coefficient [m/s]
160 const scalar alpha =
161 (2.0 + 0.6*Foam::sqrt(Re)*cbrt(Sc))*D12_/(d + ROOTVSMALL);
162
163 const scalar r = d/2;
164
165 const scalar f = F[FuelLocalId_];
166
167 const scalar rhof = composition.solids().properties()[FuelLocalId_].rho();
168
169 const scalar deltaRho0 = (nuOx_/nuFuel_)*rhof/WFuel;
170
171 // Progress variable rate
172 const scalar dfdt =
173 Aeff_*(Cb/deltaRho0)
174 /(
175 r/3/alpha
176 + sqr(r)*(1/cbrt(1-f)-1)/3/Deff
177 - (1/sqr(cbrt(1-f)))*r/k/sigma_/E_/3
178 );
179
180 // Update new progress variable
181 F[FuelLocalId_] += dfdt*dt;
182
183 // Interface radius
184 const scalar ri = r*cbrt(1-f);
185
186 // Interface radius rate
187 //const scalar dridt = -dfdt*(r/3)*pow(1-f, -2/3);
188 const scalar dridt = -dfdt*(pow3(r)/3)/sqr(ri);
189
190 // O2 flux [Kmol/sec]
191 const scalar q02 = deltaRho0*4*constant::mathematical::pi*sqr(ri)*dridt;
192
193 // Calculate the number of molar units reacted [Kmol]
194 const scalar dOmega = q02*dt;
195
196 // Heat of Reaction
197 const scalar Hc =
198 composition.solids().properties()[ProdLocalId_].Hf()
199 - composition.solids().properties()[FuelLocalId_].Hf();
200
201 //Stoichiometric mass ratio for fuel
202 const scalar sFuel = nuFuel_/(nuOx_);
203
204 //Stoichiometric mass ratio for product
205 const scalar sProd = nuProd_/(nuOx_);
206
207 // Add to carrier phase mass transfer [Kg]
208 dMassSRCarrier[O2GlobalId_] += dOmega*WO2_;
209
210 // Remove to particle mass transfer
211 dMassSolid[FuelLocalId_] -= dOmega*WFuel*sFuel;
212
213 // Add to particle product
214 dMassSolid[ProdLocalId_] += dOmega*WProd*sProd;
215
216 if (debug)
217 {
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;
227 }
228
229 // Heat of reaction [J]
230 return -dOmega*WFuel*sFuel*Hc;
231}
232
233
234template<class CloudType>
236{
237 return 1;
238}
239
240
241// ************************************************************************* //
label k
An indexed form of CGAL::Triangulation_cell_base_3<K> used to keep track of the Delaunay cells (tets)...
Definition: indexedCell.H:95
const CloudType & owner() const
Return const access to the owner cloud.
Templated base class for dsmc cloud.
Definition: DSMCCloud.H:75
Base class for heterogeneous reacting models.
Heteregeneous noncatalytic reaction MUCS approach. Reference: D. Papanastassiou and G....
virtual label nReactions() const
Number of reactions in the model.
Thermo package for (S)olids (L)iquids and (G)ases Takes reference to thermo package,...
Definition: SLGThermo.H:67
virtual scalar rho(const label speciei, const scalar p, const scalar T) const =0
Density [kg/m3].
virtual scalar W(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & coeffDict() const
Return const access to the coefficients dictionary.
Definition: subModelBase.C:131
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
basicSpecieMixture & composition
const volScalarField & T
surfaceScalarField rhof(fvc::interpolate(rho, "div(phi,rho)"))
volVectorField F(fluid.F())
constexpr scalar pi(M_PI)
dimensionedScalar exp(const dimensionedScalar &ds)
scalarField Re(const UList< complex > &cf)
Extract real component.
Definition: complexField.C:159
dimensionedSymmTensor sqr(const dimensionedVector &dv)
dimensionedScalar pow3(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensionedScalar sqrt(const dimensionedScalar &ds)
dimensionedScalar cbrt(const dimensionedScalar &ds)
prefixOSstream Pout
OSstream wrapped stdout (std::cout) with parallel prefix.
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
labelList f(nPoints)
volScalarField & nu
volScalarField & alpha
dictionary dict
const Vector< label > N(dict.get< Vector< label > >("N"))