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 = thermo_.WiValue(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  thermo_.hfiValue(speciei)*thermo_.WiValue(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  thermo_.hfiValue(speciei)*thermo_.WiValue(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 = thermo_.species()["O2"];
63  const scalar Wu = thermo_.WiValue(fuelIndex_);
64 
65  stoicRatio_ =
66  (
67  thermo_.WiValue(thermo_.defaultSpecie())
68  *specieStoichCoeffs_[thermo_.defaultSpecie()]
69  + thermo_.WiValue(O2Index)*mag(specieStoichCoeffs_[O2Index])
70  )/(Wu*mag(specieStoichCoeffs_[fuelIndex_]));
71 
72  s_ = thermo_.WiValue(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]*thermo_.WiValue(speciei);
97  }
98 
99  forAll(reaction_.rhs(), i)
100  {
101  const label speciei = reaction_.rhs()[i].index;
102  Yprod0_[speciei] = thermo_.WiValue(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<< " " << thermo_.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]*thermo_.WiValue(i)
119  /(thermo_.WiValue(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  reaction_(thermo_.species(), this->subDict("reaction")),
136  stoicRatio_(dimensionedScalar("stoicRatio", dimless, 0)),
137  s_(dimensionedScalar("s", dimless, 0)),
138  qFuel_(dimensionedScalar("qFuel", sqr(dimVelocity), 0)),
139  specieStoichCoeffs_(thermo_.species().size(), 0.0),
140  Yprod0_(thermo_.species().size(), 0.0),
141  fres_(Yprod0_.size()),
142  fuelIndex_(thermo_.species()[thermo.properties().lookup("fuel")]),
143  specieProd_(Yprod0_.size(), 1),
144  wFuel_
145  (
146  IOobject
147  (
148  this->thermo().phasePropertyName("wFuel"),
149  this->mesh().time().name(),
150  this->mesh(),
151  IOobject::NO_READ,
152  IOobject::NO_WRITE
153  ),
154  this->mesh(),
156  ),
157  semiImplicit_(readBool(this->coeffs_.lookup("semiImplicit")))
158 {
159  forAll(fres_, fresI)
160  {
161  IOobject header
162  (
163  "fres_" + thermo_.species()[fresI],
164  this->mesh().time().name(),
165  this->mesh()
166  );
167 
168  fres_.set
169  (
170  fresI,
171  new volScalarField
172  (
173  header,
174  this->mesh(),
175  dimensionedScalar("fres" + Foam::name(fresI), dimless, 0)
176  )
177  );
178  }
179 
180  calculateqFuel();
181 
183 
185 
186  if (semiImplicit_)
187  {
188  Info<< "Combustion mode: semi-implicit" << endl;
189  }
190  else
191  {
192  Info<< "Combustion mode: explicit" << endl;
193  }
194 }
195 
196 
197 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
198 
200 {}
201 
202 
203 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
204 
207 {
208  return wFuel_()*specieStoichCoeffs()[speciei];
209 }
210 
211 
214 {
215  const label specieI = thermo_.species()[Y.member()];
216 
217  volScalarField wSpecie
218  (
219  wFuel_*specieStoichCoeffs()[specieI]
220  );
221 
222  if (semiImplicit_)
223  {
224  const label fNorm = specieProd()[specieI];
225  const volScalarField fres(this->fres(specieI));
226  wSpecie /= max(fNorm*(Y - fres), scalar(1e-2));
227 
228  return -fNorm*wSpecie*fres + fNorm*fvm::Sp(wSpecie, Y);
229  }
230  else
231  {
232  return wSpecie + fvm::Sp(0.0*wSpecie, Y);
233  }
234 }
235 
236 
239 {
240  const label fuelI = fuelIndex();
241  volScalarField& YFuel =
242  const_cast<volScalarField&>(this->thermo().Y(fuelI));
243 
244  return -qFuel()*(R(YFuel) & YFuel);
245 }
246 
247 
249 {
250  if (combustionModel::read())
251  {
252  return true;
253  }
254  else
255  {
256  return false;
257  }
258 }
259 
260 
262 {
263  const label O2Index = thermo_.species()["O2"];
264  const volScalarField& YFuel = thermo_.Y()[fuelIndex_];
265  const volScalarField& YO2 = thermo_.Y()[O2Index];
266 
267  // reactants
268  forAll(reaction_.lhs(), i)
269  {
270  const label speciei = reaction_.lhs()[i].index;
271  if (speciei == fuelIndex_)
272  {
273  fres_[speciei] = max(YFuel - YO2/s_, scalar(0));
274  }
275  else if (speciei == O2Index)
276  {
277  fres_[speciei] = max(YO2 - YFuel*s_, scalar(0));
278  }
279  }
280 
281  // products
282  forAll(reaction_.rhs(), i)
283  {
284  const label speciei = reaction_.rhs()[i].index;
285  if (speciei != thermo_.defaultSpecie())
286  {
287  forAll(fres_[speciei], celli)
288  {
289  if (fres_[fuelIndex_][celli] > 0.0)
290  {
291  // rich mixture
292  fres_[speciei][celli] =
293  Yprod0_[speciei]
294  * (1.0 + YO2[celli]/s_.value() - YFuel[celli]);
295  }
296  else
297  {
298  // lean mixture
299  fres_[speciei][celli] =
300  Yprod0_[speciei]
301  * (
302  1.0
303  - YO2[celli]/s_.value()*stoicRatio_.value()
304  + YFuel[celli]*stoicRatio_.value()
305  );
306  }
307  }
308  }
309  }
310 }
311 
312 
313 // ************************************************************************* //
#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
Base class for combustion models.
const fluidMulticomponentThermo & thermo_
Reference to the thermo.
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.
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.
virtual scalar hfiValue(const label speciei) const =0
Enthalpy of formation [J/kg].
virtual const speciesTable & species() const =0
The table of species.
virtual scalar WiValue(const label speciei) const =0
Molecular weight [kg/kmol].
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 > &)
bool readBool(Istream &)
Definition: boolIO.C:66
const doubleScalar e
Definition: doubleScalar.H:106
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:257
word name(const bool)
Return a word representation of a bool.
Definition: boolIO.C:39
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
static const char nl
Definition: Ostream.H:266
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31