Skip to content

Commit 6e567d9

Browse files
dfCanteraMixture
1 parent 998dc46 commit 6e567d9

File tree

5 files changed

+711
-0
lines changed

5 files changed

+711
-0
lines changed

src/dfCanteraMixture/Make/files

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,5 @@
11
CanteraMixture.C
2+
dfSingleStepReactingMixture/dfSingleStepReactingMixture.C
23
makeThermos.C
34

45

Lines changed: 382 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,382 @@
1+
/*---------------------------------------------------------------------------*\
2+
========= |
3+
\\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4+
\\ / O peration | Website: https://openfoam.org
5+
\\ / A nd | Copyright (C) 2011-2018 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 "dfSingleStepReactingMixture.H"
27+
#include "fvMesh.H"
28+
29+
// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30+
31+
32+
void Foam::dfSingleStepReactingMixture::calculateqFuel()
33+
{
34+
// 第0条反应
35+
// const Reaction& reaction = this->operator[](0);
36+
37+
38+
// 遍历反应左边各组分
39+
// Cantera可以提供stoichCoeff
40+
// forAll(reaction.lhs(), i)
41+
// {
42+
// const label speciei = reaction.lhs()[i].index;//组分索引
43+
// const scalar stoichCoeff = reaction.lhs()[i].stoichCoeff;
44+
// specieStoichCoeffs_[speciei] = -stoichCoeff;
45+
// qFuel_.value() += this->speciesData()[speciei].hc()*stoichCoeff/Wu;
46+
// }
47+
48+
// 遍历反应右边各组分
49+
// forAll(reaction.rhs(), i)
50+
// {
51+
// const label speciei = reaction.rhs()[i].index;
52+
// const scalar stoichCoeff = reaction.rhs()[i].stoichCoeff;
53+
// specieStoichCoeffs_[speciei] = stoichCoeff;
54+
// qFuel_.value() -= this->speciesData()[speciei].hc()*stoichCoeff/Wu;
55+
// specieProd_[speciei] = -1;//反应物是1,产物是-1
56+
// }
57+
58+
const scalar Wu = this->Wi(fuelIndex_);
59+
Info << " reactant number = " << reactants_.size() << endl;
60+
forAll(reactants_, i)
61+
{
62+
const int speciei = reactants_[i];
63+
const double stoichCoeff = this->CanteraKinetics()->reactantStoichCoeff(speciei, 0);
64+
Info << "reactant stoich coeffs :" << stoichCoeff << endl;
65+
specieStoichCoeffs_[speciei] = -stoichCoeff;
66+
qFuel_.value() += this->CanteraGas()->Hf298SS(speciei)*stoichCoeff/Wu;
67+
}
68+
forAll(products_, i)
69+
{
70+
const int speciei = products_[i];
71+
const double stoichCoeff = this->CanteraKinetics()->productStoichCoeff(speciei, 0);
72+
specieStoichCoeffs_[speciei] = stoichCoeff;
73+
qFuel_.value() -= this->CanteraGas()->Hf298SS(speciei)*stoichCoeff/Wu;
74+
specieProd_[speciei] = -1;
75+
}
76+
77+
Info << "Fuel heat of combustion :" << qFuel_.value() << endl;
78+
}
79+
80+
81+
82+
void Foam::dfSingleStepReactingMixture::massAndAirStoichRatios()
83+
{
84+
const label O2Index = this->species()["O2"];
85+
// const scalar Wu = this->speciesData()[fuelIndex_].W();
86+
const scalar Wu = this->Wi(fuelIndex_);
87+
88+
stoicRatio_ =
89+
(this->Wi(inertIndex_)
90+
* specieStoichCoeffs_[inertIndex_]
91+
+ this->Wi(O2Index)
92+
* mag(specieStoichCoeffs_[O2Index]))
93+
/ (Wu*mag(specieStoichCoeffs_[fuelIndex_]));
94+
95+
s_ =
96+
(this->Wi(O2Index)
97+
* mag(specieStoichCoeffs_[O2Index]))
98+
/ (Wu*mag(specieStoichCoeffs_[fuelIndex_]));
99+
100+
101+
Info << "fuel stoich coeffs :" << specieStoichCoeffs_[fuelIndex_] << endl;
102+
103+
Info << "stoichiometric air-fuel ratio :" << stoicRatio_.value() << endl;
104+
105+
Info << "stoichiometric oxygen-fuel ratio :" << s_.value() << endl;
106+
}
107+
108+
109+
110+
void Foam::dfSingleStepReactingMixture::calculateMaxProducts()
111+
{
112+
// const Reaction& reaction = this->operator[](0);
113+
114+
scalar Wm = 0.0;
115+
scalar totalMol = 0.0;
116+
// forAll(reaction.rhs(), i)
117+
// {
118+
// label speciei = reaction.rhs()[i].index;
119+
// totalMol += mag(specieStoichCoeffs_[speciei]);
120+
// }
121+
std::string products = this->CanteraKinetics()->productString(0);
122+
forAll(products_, i)
123+
{
124+
const int speciei = products_[i];
125+
totalMol += mag(specieStoichCoeffs_[speciei]);
126+
}
127+
128+
// scalarList Xi(reaction.rhs().size());
129+
scalarList Xi(products_.size());
130+
131+
// forAll(reaction.rhs(), i)
132+
// {
133+
// const label speciei = reaction.rhs()[i].index;
134+
// Xi[i] = mag(specieStoichCoeffs_[speciei])/totalMol;
135+
136+
// Wm += Xi[i]*this->speciesData()[speciei].W();
137+
// }
138+
139+
forAll(products_, i)
140+
{
141+
const int speciei = products_[i];
142+
Xi[i] = mag(specieStoichCoeffs_[speciei])/totalMol;
143+
144+
Wm += Xi[i]*this->Wi(speciei);
145+
}
146+
147+
// forAll(reaction.rhs(), i)
148+
// {
149+
// const label speciei = reaction.rhs()[i].index;
150+
// Yprod0_[speciei] = this->speciesData()[speciei].W()/Wm*Xi[i];
151+
// }
152+
153+
forAll(products_, i)
154+
{
155+
const int speciei = products_[i];
156+
Yprod0_[speciei] = this->Wi(speciei)/Wm*Xi[i];
157+
}
158+
159+
Info << "Maximum products mass concentrations:" << nl;
160+
forAll(Yprod0_, i)
161+
{
162+
if (Yprod0_[i] > 0)
163+
{
164+
Info<< " " << this->species()[i] << ": " << Yprod0_[i] << nl;
165+
}
166+
}
167+
168+
// Normalize the stoichiometric coeff to mass
169+
forAll(specieStoichCoeffs_, i)
170+
{
171+
specieStoichCoeffs_[i] =
172+
specieStoichCoeffs_[i]
173+
* this->Wi(i)
174+
/ (this->Wi(fuelIndex_)
175+
* mag(specieStoichCoeffs_[fuelIndex_]));
176+
}
177+
}
178+
179+
180+
181+
void Foam::dfSingleStepReactingMixture::fresCorrect()
182+
{
183+
// const Reaction& reaction = this->operator[](0);
184+
185+
label O2Index = this->species()["O2"];
186+
const volScalarField& YFuel = this->Y()[fuelIndex_];
187+
const volScalarField& YO2 = this->Y()[O2Index];
188+
189+
// // reactants
190+
// forAll(reaction.lhs(), i)
191+
// {
192+
// const label speciei = reaction.lhs()[i].index;
193+
// if (speciei == fuelIndex_)
194+
// {
195+
// fres_[speciei] = max(YFuel - YO2/s_, scalar(0));
196+
// }
197+
// else if (speciei == O2Index)
198+
// {
199+
// fres_[speciei] = max(YO2 - YFuel*s_, scalar(0));
200+
// }
201+
// }
202+
203+
forAll(reactants_, i)
204+
{
205+
const int speciei = reactants_[i];
206+
if (speciei == fuelIndex_)
207+
{
208+
fres_[speciei] = max(YFuel - YO2/s_, scalar(0));
209+
}
210+
else if (speciei == O2Index)
211+
{
212+
fres_[speciei] = max(YO2 - YFuel*s_, scalar(0));
213+
}
214+
}
215+
216+
217+
// // products
218+
// forAll(reaction.rhs(), i)
219+
// {
220+
// const label speciei = reaction.rhs()[i].index;
221+
// if (speciei != inertIndex_)
222+
// {
223+
// forAll(fres_[speciei], celli)
224+
// {
225+
// if (fres_[fuelIndex_][celli] > 0.0)
226+
// {
227+
// // rich mixture
228+
// fres_[speciei][celli] =
229+
// Yprod0_[speciei]
230+
// * (1.0 + YO2[celli]/s_.value() - YFuel[celli]);
231+
// }
232+
// else
233+
// {
234+
// // lean mixture
235+
// fres_[speciei][celli] =
236+
// Yprod0_[speciei]
237+
// * (
238+
// 1.0
239+
// - YO2[celli]/s_.value()*stoicRatio_.value()
240+
// + YFuel[celli]*stoicRatio_.value()
241+
// );
242+
// }
243+
// }
244+
// }
245+
// }
246+
forAll(products_, i)
247+
{
248+
const int speciei = products_[i];
249+
if (speciei != inertIndex_)
250+
{
251+
forAll(fres_[speciei], celli)
252+
{
253+
if (fres_[fuelIndex_][celli] > 0.0)
254+
{
255+
// rich mixture
256+
fres_[speciei][celli] =
257+
Yprod0_[speciei]
258+
* (1.0 + YO2[celli]/s_.value() - YFuel[celli]);
259+
}
260+
else
261+
{
262+
// lean mixture
263+
fres_[speciei][celli] =
264+
Yprod0_[speciei]
265+
* (
266+
1.0
267+
- YO2[celli]/s_.value()*stoicRatio_.value()
268+
+ YFuel[celli]*stoicRatio_.value()
269+
);
270+
}
271+
}
272+
}
273+
}
274+
}
275+
276+
277+
// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
278+
279+
280+
Foam::dfSingleStepReactingMixture::dfSingleStepReactingMixture
281+
(
282+
const dictionary& thermoDict,
283+
const fvMesh& mesh,
284+
const word& phaseName
285+
)
286+
:
287+
CanteraMixture(thermoDict, mesh, phaseName),
288+
stoicRatio_(dimensionedScalar("stoicRatio", dimless, 0)),
289+
s_(dimensionedScalar("s", dimless, 0)),
290+
qFuel_(dimensionedScalar("qFuel", sqr(dimVelocity), 0)),
291+
specieStoichCoeffs_(this->species().size(), 0.0),
292+
Yprod0_(this->species().size(), 0.0),
293+
fres_(Yprod0_.size()),
294+
inertIndex_(this->species()[thermoDict.lookup("inertSpecie")]),
295+
fuelIndex_(this->species()[thermoDict.lookup("fuel")]),
296+
specieProd_(Yprod0_.size(), 1)
297+
{
298+
if (this->nReactions() == 1)
299+
{
300+
forAll(fres_, fresI)
301+
{
302+
IOobject header
303+
(
304+
"fres_" + this->species()[fresI],
305+
mesh.time().timeName(),
306+
mesh,
307+
IOobject::NO_READ,
308+
IOobject::NO_WRITE
309+
);
310+
311+
fres_.set
312+
(
313+
fresI,
314+
new volScalarField
315+
(
316+
header,
317+
mesh,
318+
dimensionedScalar("fres" + name(fresI), dimless, 0)
319+
)
320+
);
321+
}
322+
323+
std::string reactants = this->CanteraKinetics()->reactantString(0);
324+
std::string products = this->CanteraKinetics()->productString(0);
325+
Info << "reactants: " << reactants << endl;
326+
for(int i=0; i<this->nSpecies(); ++i)
327+
{
328+
word name = this->CanteraGas()->speciesName(i);
329+
if(reactants.find(name) != std::string::npos)
330+
{
331+
reactants_.append(i);
332+
}
333+
else if(products.find(name) != std::string::npos)
334+
{
335+
products_.append(i);
336+
}
337+
}
338+
Info << "reactants: " << reactants_ << endl;
339+
340+
calculateqFuel();
341+
342+
massAndAirStoichRatios();
343+
344+
calculateMaxProducts();
345+
}
346+
else
347+
{
348+
FatalErrorInFunction
349+
<< "Only one reaction required for single step reaction"
350+
<< exit(FatalError);
351+
}
352+
353+
// std::string reactants = this->CanteraKinetics()->reactantString(0);
354+
// std::string products = this->CanteraKinetics()->productString(0);
355+
// Info << "reactants: " << reactants << endl;
356+
// for(int i=0; i<this->nSpecies(); ++i)
357+
// {
358+
// word name = this->CanteraGas()->speciesName(i);
359+
// if(reactants.find(name) != std::string::npos)
360+
// {
361+
// reactants_.append(i);
362+
// }
363+
// else if(products.find(name) != std::string::npos)
364+
// {
365+
// products_.append(i);
366+
// }
367+
// }
368+
// Info << "reactants: " << reactants_ << endl;
369+
}
370+
371+
372+
// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
373+
374+
375+
void Foam::dfSingleStepReactingMixture::read
376+
(
377+
const dictionary& thermoDict
378+
)
379+
{}
380+
381+
382+
// ************************************************************************* //

0 commit comments

Comments
 (0)