singleStepReactingMixture.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-------------------------------------------------------------------------------
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
29#include "fvMesh.H"
30
31// * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32
33template<class ThermoType>
35{
36 const Reaction<ThermoType>& reaction = this->operator[](0);
37 const scalar Wu = this->speciesData()[fuelIndex_].W();
38
39 forAll(reaction.lhs(), i)
40 {
41 const label speciei = reaction.lhs()[i].index;
42 const scalar stoichCoeff = reaction.lhs()[i].stoichCoeff;
43 specieStoichCoeffs_[speciei] = -stoichCoeff;
44 qFuel_.value() += this->speciesData()[speciei].hc()*stoichCoeff/Wu;
45 }
46
47 forAll(reaction.rhs(), i)
48 {
49 const label speciei = reaction.rhs()[i].index;
50 const scalar stoichCoeff = reaction.rhs()[i].stoichCoeff;
51 specieStoichCoeffs_[speciei] = stoichCoeff;
52 qFuel_.value() -= this->speciesData()[speciei].hc()*stoichCoeff/Wu;
53 specieProd_[speciei] = -1;
54 }
55
56 Info<< "Fuel heat of combustion :" << qFuel_.value() << endl;
57}
58
59
60template<class ThermoType>
62{
63 const label O2Index = this->species().find("O2");
64 const scalar Wu = this->speciesData()[fuelIndex_].W();
65
66 stoicRatio_ =
67 (this->speciesData()[inertIndex_].W()
68 * specieStoichCoeffs_[inertIndex_]
69 + this->speciesData()[O2Index].W()
70 * mag(specieStoichCoeffs_[O2Index]))
71 / (Wu*mag(specieStoichCoeffs_[fuelIndex_]));
72
73 s_ =
74 (this->speciesData()[O2Index].W()
75 * mag(specieStoichCoeffs_[O2Index]))
76 / (Wu*mag(specieStoichCoeffs_[fuelIndex_]));
77
78 Info<< "stoichiometric air-fuel ratio :" << stoicRatio_.value() << endl;
79
80 Info<< "stoichiometric oxygen-fuel ratio :" << s_.value() << endl;
81}
82
83
84template<class ThermoType>
86{
87 const Reaction<ThermoType>& reaction = this->operator[](0);
88
89 scalar Wm = 0.0;
90 scalar totalMol = 0.0;
91 forAll(reaction.rhs(), i)
92 {
93 label speciei = reaction.rhs()[i].index;
94 totalMol += mag(specieStoichCoeffs_[speciei]);
95 }
96
97 scalarList Xi(reaction.rhs().size());
98
99 forAll(reaction.rhs(), i)
100 {
101 const label speciei = reaction.rhs()[i].index;
102 Xi[i] = mag(specieStoichCoeffs_[speciei])/totalMol;
103
104 Wm += Xi[i]*this->speciesData()[speciei].W();
105 }
106
107 forAll(reaction.rhs(), i)
108 {
109 const label speciei = reaction.rhs()[i].index;
110 Yprod0_[speciei] = this->speciesData()[speciei].W()/Wm*Xi[i];
111 }
112
113 Info<< "Maximum products mass concentrations:" << nl;
114 forAll(Yprod0_, i)
115 {
116 if (Yprod0_[i] > 0)
117 {
118 Info<< " " << this->species()[i] << ": " << Yprod0_[i] << nl;
119 }
120 }
121
122 // Normalize the stoichiometric coeff to mass
123 forAll(specieStoichCoeffs_, i)
124 {
125 specieStoichCoeffs_[i] =
126 specieStoichCoeffs_[i]
127 * this->speciesData()[i].W()
128 / (this->speciesData()[fuelIndex_].W()
129 * mag(specieStoichCoeffs_[fuelIndex_]));
130 }
131}
132
133
134template<class ThermoType>
136{
137 const Reaction<ThermoType>& reaction = this->operator[](0);
138
139 const label O2Index = this->species().find("O2");
140 const volScalarField& YFuel = this->Y()[fuelIndex_];
141 const volScalarField& YO2 = this->Y()[O2Index];
142
143 // reactants
144 forAll(reaction.lhs(), i)
145 {
146 const label speciei = reaction.lhs()[i].index;
147 if (speciei == fuelIndex_)
148 {
149 fres_[speciei] = max(YFuel - YO2/s_, scalar(0));
150 }
151 else if (speciei == O2Index)
152 {
153 fres_[speciei] = max(YO2 - YFuel*s_, scalar(0));
154 }
155 }
156
157
158 // products
159 forAll(reaction.rhs(), i)
160 {
161 const label speciei = reaction.rhs()[i].index;
162 if (speciei != inertIndex_)
163 {
164 forAll(fres_[speciei], celli)
165 {
166 if (fres_[fuelIndex_][celli] > 0.0)
167 {
168 // rich mixture
169 fres_[speciei][celli] =
170 Yprod0_[speciei]
171 * (1.0 + YO2[celli]/s_.value() - YFuel[celli]);
172 }
173 else
174 {
175 // lean mixture
176 fres_[speciei][celli] =
177 Yprod0_[speciei]
178 * (
179 1.0
180 - YO2[celli]/s_.value()*stoicRatio_.value()
181 + YFuel[celli]*stoicRatio_.value()
182 );
183 }
184 }
185 }
186 }
187}
188
189
190// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
191
192template<class ThermoType>
194(
195 const dictionary& thermoDict,
196 const fvMesh& mesh,
197 const word& phaseName
198)
199:
200 reactingMixture<ThermoType>(thermoDict, mesh, phaseName),
201 stoicRatio_(dimensionedScalar("stoicRatio", dimless, Zero)),
202 s_(dimensionedScalar("s", dimless, Zero)),
203 qFuel_(dimensionedScalar("qFuel", sqr(dimVelocity), Zero)),
204 specieStoichCoeffs_(this->species_.size(), Zero),
205 Yprod0_(this->species_.size(), Zero),
206 fres_(Yprod0_.size()),
207 inertIndex_(this->species().find(thermoDict.get<word>("inertSpecie"))),
208 fuelIndex_(this->species().find(thermoDict.get<word>("fuel"))),
209 specieProd_(Yprod0_.size(), 1)
210{
211 if (this->size() == 1)
212 {
213 forAll(fres_, fresI)
214 {
215 IOobject header
216 (
217 "fres_" + this->species()[fresI],
218 mesh.time().timeName(),
219 mesh,
222 );
223
224 fres_.set
225 (
226 fresI,
228 (
229 header,
230 mesh,
232 )
233 );
234 }
235
237
239
241 }
242 else
243 {
245 << "Only one reaction required for single step reaction"
246 << exit(FatalError);
247 }
248}
249
250
251// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
252
253template<class ThermoType>
255(
257)
258{}
259
260
261// ************************************************************************* //
Y[inertIndex] max(0.0)
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:170
virtual bool read()
Re-read model coefficients if they have changed.
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:73
static word timeName(const scalar t, const int precision=precision_)
Definition: Time.C:780
label size() const noexcept
The number of elements in the list.
Definition: UPtrListI.H:106
const speciesTable & species() const
Return the table of species.
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:91
const Time & time() const
Return the top-level database.
Definition: fvMesh.H:290
Foam::reactingMixture.
Single step reacting mixture.
void fresCorrect()
Calculates the residual for all components.
void calculateMaxProducts()
Calculate maximum products at stoichiometric mixture.
PtrList< volScalarField > fres_
List of components residual.
void massAndAirStoichRatios()
Calculate air/fuel and oxygen/fuel ratio.
A class for handling words, derived from Foam::string.
Definition: word.H:68
PtrList< volScalarField > & Y
const dictionary & thermoDict
Definition: EEqn.H:16
dynamicFvMesh & mesh
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
CombustionModel< rhoReactionThermo > & reaction
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
const dimensionSet dimless
Dimensionless.
dimensionedSymmTensor sqr(const dimensionedVector &dv)
const dimensionSet dimVelocity
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
constexpr char nl
The newline '\n' character (0x0a)
Definition: Ostream.H:53
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333