singleStepReactingMixture.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) 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 
27 #include "fvMesh.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 {
34  const Reaction<ThermoType>& reaction = this->operator[](0);
35  const scalar Wu = this->speciesData()[fuelIndex_].W();
36 
37  forAll(reaction.lhs(), i)
38  {
39  const label speciei = reaction.lhs()[i].index;
40  const scalar stoichCoeff = reaction.lhs()[i].stoichCoeff;
41  specieStoichCoeffs_[speciei] = -stoichCoeff;
42  qFuel_.value() += this->speciesData()[speciei].hc()*stoichCoeff/Wu;
43  }
44 
45  forAll(reaction.rhs(), i)
46  {
47  const label speciei = reaction.rhs()[i].index;
48  const scalar stoichCoeff = reaction.rhs()[i].stoichCoeff;
49  specieStoichCoeffs_[speciei] = stoichCoeff;
50  qFuel_.value() -= this->speciesData()[speciei].hc()*stoichCoeff/Wu;
51  specieProd_[speciei] = -1;
52  }
53 
54  Info << "Fuel heat of combustion :" << qFuel_.value() << endl;
55 }
56 
57 
58 template<class ThermoType>
60 {
61  const label O2Index = this->species()["O2"];
62  const scalar Wu = this->speciesData()[fuelIndex_].W();
63 
64  stoicRatio_ =
65  (this->speciesData()[inertIndex_].W()
66  * specieStoichCoeffs_[inertIndex_]
67  + this->speciesData()[O2Index].W()
68  * mag(specieStoichCoeffs_[O2Index]))
69  / (Wu*mag(specieStoichCoeffs_[fuelIndex_]));
70 
71  s_ =
72  (this->speciesData()[O2Index].W()
73  * mag(specieStoichCoeffs_[O2Index]))
74  / (Wu*mag(specieStoichCoeffs_[fuelIndex_]));
75 
76  Info << "stoichiometric air-fuel ratio :" << stoicRatio_.value() << endl;
77 
78  Info << "stoichiometric oxygen-fuel ratio :" << s_.value() << endl;
79 }
80 
81 
82 template<class ThermoType>
84 {
85  const Reaction<ThermoType>& reaction = this->operator[](0);
86 
87  scalar Wm = 0.0;
88  scalar totalMol = 0.0;
89  forAll(reaction.rhs(), i)
90  {
91  label speciei = reaction.rhs()[i].index;
92  totalMol += mag(specieStoichCoeffs_[speciei]);
93  }
94 
95  scalarList Xi(reaction.rhs().size());
96 
97  forAll(reaction.rhs(), i)
98  {
99  const label speciei = reaction.rhs()[i].index;
100  Xi[i] = mag(specieStoichCoeffs_[speciei])/totalMol;
101 
102  Wm += Xi[i]*this->speciesData()[speciei].W();
103  }
104 
105  forAll(reaction.rhs(), i)
106  {
107  const label speciei = reaction.rhs()[i].index;
108  Yprod0_[speciei] = this->speciesData()[speciei].W()/Wm*Xi[i];
109  }
110 
111  Info << "Maximum products mass concentrations:" << nl;
112  forAll(Yprod0_, i)
113  {
114  if (Yprod0_[i] > 0)
115  {
116  Info<< " " << this->species()[i] << ": " << Yprod0_[i] << nl;
117  }
118  }
119 
120  // Normalize the stoichiometric coeff to mass
121  forAll(specieStoichCoeffs_, i)
122  {
123  specieStoichCoeffs_[i] =
124  specieStoichCoeffs_[i]
125  * this->speciesData()[i].W()
126  / (this->speciesData()[fuelIndex_].W()
127  * mag(specieStoichCoeffs_[fuelIndex_]));
128  }
129 }
130 
131 
132 template<class ThermoType>
134 {
135  const Reaction<ThermoType>& reaction = this->operator[](0);
136 
137  label O2Index = this->species()["O2"];
138  const volScalarField& YFuel = this->Y()[fuelIndex_];
139  const volScalarField& YO2 = this->Y()[O2Index];
140 
141  // reactants
142  forAll(reaction.lhs(), i)
143  {
144  const label speciei = reaction.lhs()[i].index;
145  if (speciei == fuelIndex_)
146  {
147  fres_[speciei] = max(YFuel - YO2/s_, scalar(0));
148  }
149  else if (speciei == O2Index)
150  {
151  fres_[speciei] = max(YO2 - YFuel*s_, scalar(0));
152  }
153  }
154 
155 
156  // products
157  forAll(reaction.rhs(), i)
158  {
159  const label speciei = reaction.rhs()[i].index;
160  if (speciei != inertIndex_)
161  {
162  forAll(fres_[speciei], celli)
163  {
164  if (fres_[fuelIndex_][celli] > 0.0)
165  {
166  // rich mixture
167  fres_[speciei][celli] =
168  Yprod0_[speciei]
169  * (1.0 + YO2[celli]/s_.value() - YFuel[celli]);
170  }
171  else
172  {
173  // lean mixture
174  fres_[speciei][celli] =
175  Yprod0_[speciei]
176  * (
177  1.0
178  - YO2[celli]/s_.value()*stoicRatio_.value()
179  + YFuel[celli]*stoicRatio_.value()
180  );
181  }
182  }
183  }
184  }
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
189 
190 template<class ThermoType>
192 (
193  const dictionary& thermoDict,
194  const fvMesh& mesh,
195  const word& phaseName
196 )
197 :
198  reactingMixture<ThermoType>(thermoDict, mesh, phaseName),
199  stoicRatio_(dimensionedScalar("stoicRatio", dimless, 0)),
200  s_(dimensionedScalar("s", dimless, 0)),
201  qFuel_(dimensionedScalar("qFuel", sqr(dimVelocity), 0)),
202  specieStoichCoeffs_(this->species_.size(), 0.0),
203  Yprod0_(this->species_.size(), 0.0),
204  fres_(Yprod0_.size()),
205  inertIndex_(this->species()[thermoDict.lookup("inertSpecie")]),
206  fuelIndex_(this->species()[thermoDict.lookup("fuel")]),
207  specieProd_(Yprod0_.size(), 1)
208 {
209  if (this->size() == 1)
210  {
211  forAll(fres_, fresI)
212  {
213  IOobject header
214  (
215  "fres_" + this->species()[fresI],
216  mesh.time().timeName(),
217  mesh,
218  IOobject::NO_READ,
219  IOobject::NO_WRITE
220  );
221 
222  fres_.set
223  (
224  fresI,
225  new volScalarField
226  (
227  header,
228  mesh,
229  dimensionedScalar("fres" + name(fresI), dimless, 0)
230  )
231  );
232  }
233 
234  calculateqFuel();
235 
236  massAndAirStoichRatios();
237 
238  calculateMaxProducts();
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 
253 template<class ThermoType>
255 (
256  const dictionary& thermoDict
257 )
258 {}
259 
260 
261 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Foam::reactingMixture.
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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
A list of keyword definitions, which are a keyword followed by any number of values (e...
Definition: dictionary.H:158
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:319
tUEqn clear()
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
const List< specieCoeffs > & lhs() const
Return the components of the left hand side.
Definition: ReactionI.H:57
void massAndAirStoichRatios()
Calculate air/fuel and oxygen/fuel ratio.
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:56
dynamicFvMesh & mesh
A class for handling words, derived from string.
Definition: word.H:59
CombustionModel< rhoReactionThermo > & reaction
static const char nl
Definition: Ostream.H:265
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
void calculateMaxProducts()
Calculate maximum products at stoichiometric mixture.
void read(const dictionary &)
Read dictionary.
const dimensionSet dimless(0, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:47
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
singleStepReactingMixture(const dictionary &, const fvMesh &, const word &)
Construct from dictionary, mesh and phase name.
messageStream Info
dimensioned< scalar > mag(const dimensioned< Type > &)
An auto-pointer similar to the STL auto_ptr but with automatic casting to a reference to the type and...
Definition: PtrList.H:52
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
Definition: ReactionI.H:64
void fresCorrect()
Calculates the residual for all components.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
const dimensionSet dimVelocity