singleStepCombustion.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-2023 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 "fvSolution.H"
27 #include "singleStepCombustion.H"
28 #include "fvmSup.H"
29 
30 
31 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
32 
34 {
35  const scalar Wu = mixture_.Wi(fuelIndex_);
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() +=
43  mixture_.Hf(speciei)*mixture_.Wi(speciei)*stoichCoeff/Wu;
44  }
45 
46  forAll(reaction_.rhs(), i)
47  {
48  const label speciei = reaction_.rhs()[i].index;
49  const scalar stoichCoeff = reaction_.rhs()[i].stoichCoeff;
50  specieStoichCoeffs_[speciei] = stoichCoeff;
51  qFuel_.value() -=
52  mixture_.Hf(speciei)*mixture_.Wi(speciei)*stoichCoeff/Wu;
53  specieProd_[speciei] = -1;
54  }
55 
56  Info << "Fuel heat of combustion: " << qFuel_.value() << endl;
57 }
58 
59 
61 {
62  const label O2Index = mixture_.species()["O2"];
63  const scalar Wu = mixture_.Wi(fuelIndex_);
64 
65  stoicRatio_ =
66  (
67  mixture_.Wi(mixture_.defaultSpecie())
68  *specieStoichCoeffs_[mixture_.defaultSpecie()]
69  + mixture_.Wi(O2Index)*mag(specieStoichCoeffs_[O2Index])
70  )/(Wu*mag(specieStoichCoeffs_[fuelIndex_]));
71 
72  s_ = mixture_.Wi(O2Index)*mag(specieStoichCoeffs_[O2Index])
73  /(Wu*mag(specieStoichCoeffs_[fuelIndex_]));
74 
75  Info << "stoichiometric air-fuel ratio: " << stoicRatio_.value() << endl;
76  Info << "stoichiometric oxygen-fuel ratio: " << s_.value() << endl;
77 }
78 
79 
81 {
82  scalar Wm = 0.0;
83  scalar totalMol = 0.0;
84  forAll(reaction_.rhs(), i)
85  {
86  label speciei = reaction_.rhs()[i].index;
87  totalMol += mag(specieStoichCoeffs_[speciei]);
88  }
89 
90  scalarList Xi(reaction_.rhs().size());
91 
92  forAll(reaction_.rhs(), i)
93  {
94  const label speciei = reaction_.rhs()[i].index;
95  Xi[i] = mag(specieStoichCoeffs_[speciei])/totalMol;
96  Wm += Xi[i]*mixture_.Wi(speciei);
97  }
98 
99  forAll(reaction_.rhs(), i)
100  {
101  const label speciei = reaction_.rhs()[i].index;
102  Yprod0_[speciei] = mixture_.Wi(speciei)/Wm*Xi[i];
103  }
104 
105  Info << "Maximum products mass concentrations: " << nl;
106  forAll(Yprod0_, i)
107  {
108  if (Yprod0_[i] > 0)
109  {
110  Info<< " " << mixture_.species()[i] << ": " << Yprod0_[i] << nl;
111  }
112  }
113 
114  // Normalise the stoichiometric coeff to mass
115  forAll(specieStoichCoeffs_, i)
116  {
117  specieStoichCoeffs_[i] =
118  specieStoichCoeffs_[i]*mixture_.Wi(i)
119  /(mixture_.Wi(fuelIndex_)*mag(specieStoichCoeffs_[fuelIndex_]));
120  }
121 }
122 
123 
124 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
125 
127 (
128  const word& modelType,
131  const word& combustionProperties
132 )
133 :
134  combustionModel(modelType, thermo, turb, combustionProperties),
135  mixture_(dynamic_cast<const basicSpecieMixture&>(this->thermo())),
136  reaction_(mixture_.species(), this->subDict("reaction")),
137  stoicRatio_(dimensionedScalar("stoicRatio", dimless, 0)),
138  s_(dimensionedScalar("s", dimless, 0)),
139  qFuel_(dimensionedScalar("qFuel", sqr(dimVelocity), 0)),
140  specieStoichCoeffs_(mixture_.species().size(), 0.0),
141  Yprod0_(mixture_.species().size(), 0.0),
142  fres_(Yprod0_.size()),
143  fuelIndex_(mixture_.species()[thermo.properties().lookup("fuel")]),
144  specieProd_(Yprod0_.size(), 1),
145  wFuel_
146  (
147  IOobject
148  (
149  this->thermo().phasePropertyName("wFuel"),
150  this->mesh().time().name(),
151  this->mesh(),
152  IOobject::NO_READ,
153  IOobject::NO_WRITE
154  ),
155  this->mesh(),
157  ),
158  semiImplicit_(readBool(this->coeffs_.lookup("semiImplicit")))
159 {
160  forAll(fres_, fresI)
161  {
162  IOobject header
163  (
164  "fres_" + mixture_.species()[fresI],
165  this->mesh().time().name(),
166  this->mesh()
167  );
168 
169  fres_.set
170  (
171  fresI,
172  new volScalarField
173  (
174  header,
175  this->mesh(),
176  dimensionedScalar("fres" + Foam::name(fresI), dimless, 0)
177  )
178  );
179  }
180 
181  calculateqFuel();
182 
184 
186 
187  if (semiImplicit_)
188  {
189  Info<< "Combustion mode: semi-implicit" << endl;
190  }
191  else
192  {
193  Info<< "Combustion mode: explicit" << endl;
194  }
195 }
196 
197 
198 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
199 
201 {}
202 
203 
204 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
205 
208 {
209  return wFuel_()*specieStoichCoeffs()[speciei];
210 }
211 
212 
215 {
216  const label specieI = mixture_.species()[Y.member()];
217 
218  volScalarField wSpecie
219  (
220  wFuel_*specieStoichCoeffs()[specieI]
221  );
222 
223  if (semiImplicit_)
224  {
225  const label fNorm = specieProd()[specieI];
226  const volScalarField fres(this->fres(specieI));
227  wSpecie /= max(fNorm*(Y - fres), scalar(1e-2));
228 
229  return -fNorm*wSpecie*fres + fNorm*fvm::Sp(wSpecie, Y);
230  }
231  else
232  {
233  return wSpecie + fvm::Sp(0.0*wSpecie, Y);
234  }
235 }
236 
237 
240 {
241  const label fuelI = fuelIndex();
242  volScalarField& YFuel =
243  const_cast<volScalarField&>(this->thermo().composition().Y(fuelI));
244 
245  return -qFuel()*(R(YFuel) & YFuel);
246 }
247 
248 
250 {
251  if (combustionModel::read())
252  {
253  return true;
254  }
255  else
256  {
257  return false;
258  }
259 }
260 
261 
263 {
264  const label O2Index = mixture_.species()["O2"];
265  const volScalarField& YFuel = mixture_.Y()[fuelIndex_];
266  const volScalarField& YO2 = mixture_.Y()[O2Index];
267 
268  // reactants
269  forAll(reaction_.lhs(), i)
270  {
271  const label speciei = reaction_.lhs()[i].index;
272  if (speciei == fuelIndex_)
273  {
274  fres_[speciei] = max(YFuel - YO2/s_, scalar(0));
275  }
276  else if (speciei == O2Index)
277  {
278  fres_[speciei] = max(YO2 - YFuel*s_, scalar(0));
279  }
280  }
281 
282  // products
283  forAll(reaction_.rhs(), i)
284  {
285  const label speciei = reaction_.rhs()[i].index;
286  if (speciei != mixture_.defaultSpecie())
287  {
288  forAll(fres_[speciei], celli)
289  {
290  if (fres_[fuelIndex_][celli] > 0.0)
291  {
292  // rich mixture
293  fres_[speciei][celli] =
294  Yprod0_[speciei]
295  * (1.0 + YO2[celli]/s_.value() - YFuel[celli]);
296  }
297  else
298  {
299  // lean mixture
300  fres_[speciei][celli] =
301  Yprod0_[speciei]
302  * (
303  1.0
304  - YO2[celli]/s_.value()*stoicRatio_.value()
305  + YFuel[celli]*stoicRatio_.value()
306  );
307  }
308  }
309  }
310  }
311 }
312 
313 
314 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Generic GeometricField class.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
Specialisation of basicMixture for a mixture consisting of a number for molecular species.
virtual scalar Wi(const label speciei) const =0
Molecular weight of the given specie [kg/kmol].
const speciesTable & species() const
Return the table of species.
virtual scalar Hf(const label speciei) const =0
Enthalpy of formation [J/kg].
Base class for combustion models.
const fvMesh & mesh() const
Return const access to the mesh database.
virtual bool read()
Update properties from given dictionary.
void fresCorrect()
Calculates the residual for all components.
void calculateMaxProducts()
Calculate maximum products at stoichiometric mixture.
reaction reaction_
The single-step reaction.
const basicSpecieMixture & mixture_
Reference to the mixture.
List< int > specieProd_
List to indicate if specie is produced/consumed.
PtrList< volScalarField > fres_
List of components residual.
singleStepCombustion(const word &modelType, const fluidMulticomponentThermo &thermo, const compressibleMomentumTransportModel &turb, const word &combustionProperties)
Construct from components.
void massAndAirStoichRatios()
Calculate air/fuel and oxygen/fuel ratio.
scalarList specieStoichCoeffs_
Stoichiometric coefficient for the reaction.
virtual tmp< volScalarField::Internal > R(const label speciei) const
Specie consumption rate field.
bool semiImplicit_
Semi-implicit (true) or explicit (false) treatment.
virtual tmp< volScalarField > Qdot() const
Heat release rate [kg/m/s^3].
dimensionedScalar qFuel_
Heat of combustion [J/Kg].
virtual bool read()
Update properties from given dictionary.
Base class for single-phase compressible turbulence models.
const Type & value() const
Return const reference to value.
Base-class for multi-component fluid thermodynamic properties.
const List< specieCoeffs > & lhs() const
Return the components of the left hand side.
Definition: reactionI.H:36
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
Definition: reactionI.H:42
A class for managing temporary objects.
Definition: tmp.H:55
A class for handling words, derived from string.
Definition: word.H:62
Calculate the matrix for implicit and explicit sources.
tmp< fvMatrix< Type > > Sp(const volScalarField::Internal &, const VolField< Type > &)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
bool readBool(Istream &)
Definition: boolIO.C:60
const doubleScalar e
Definition: doubleScalar.H:105
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
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const dimensionSet dimless
messageStream Info
const dimensionSet dimTime
const dimensionSet dimVolume
dimensioned< scalar > mag(const dimensioned< Type > &)
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
const dimensionSet dimMass
const dimensionSet dimVelocity
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
static const char nl
Definition: Ostream.H:260
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31