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-2024 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 "singleStepCombustion.H"
27 #include "fvmSup.H"
28 
29 // * * * * * * * * * * * * Protected Member Functions * * * * * * * * * * * //
30 
32 {
33  const scalar Wu = thermo_.WiValue(fuelIndex_);
34 
35  forAll(reaction_.lhs(), i)
36  {
37  const label speciei = reaction_.lhs()[i].index;
38  const scalar stoichCoeff = reaction_.lhs()[i].stoichCoeff;
39  specieStoichCoeffs_[speciei] = -stoichCoeff;
40  qFuel_.value() +=
41  thermo_.hfiValue(speciei)*thermo_.WiValue(speciei)*stoichCoeff/Wu;
42  }
43 
44  forAll(reaction_.rhs(), i)
45  {
46  const label speciei = reaction_.rhs()[i].index;
47  const scalar stoichCoeff = reaction_.rhs()[i].stoichCoeff;
48  specieStoichCoeffs_[speciei] = stoichCoeff;
49  qFuel_.value() -=
50  thermo_.hfiValue(speciei)*thermo_.WiValue(speciei)*stoichCoeff/Wu;
51  specieProd_[speciei] = -1;
52  }
53 
54  Info << "Fuel heat of combustion: " << qFuel_.value() << endl;
55 }
56 
57 
59 {
60  const label O2Index = thermo_.species()["O2"];
61  const scalar Wu = thermo_.WiValue(fuelIndex_);
62 
63  stoicRatio_ =
64  (
65  thermo_.WiValue(thermo_.defaultSpecie())
66  *specieStoichCoeffs_[thermo_.defaultSpecie()]
67  + thermo_.WiValue(O2Index)*mag(specieStoichCoeffs_[O2Index])
68  )/(Wu*mag(specieStoichCoeffs_[fuelIndex_]));
69 
70  s_ = thermo_.WiValue(O2Index)*mag(specieStoichCoeffs_[O2Index])
71  /(Wu*mag(specieStoichCoeffs_[fuelIndex_]));
72 
73  Info << "stoichiometric air-fuel ratio: " << stoicRatio_.value() << endl;
74  Info << "stoichiometric oxygen-fuel ratio: " << s_.value() << endl;
75 }
76 
77 
79 {
80  scalar Wm = 0.0;
81  scalar totalMol = 0.0;
82  forAll(reaction_.rhs(), i)
83  {
84  label speciei = reaction_.rhs()[i].index;
85  totalMol += mag(specieStoichCoeffs_[speciei]);
86  }
87 
88  scalarList Xi(reaction_.rhs().size());
89 
90  forAll(reaction_.rhs(), i)
91  {
92  const label speciei = reaction_.rhs()[i].index;
93  Xi[i] = mag(specieStoichCoeffs_[speciei])/totalMol;
94  Wm += Xi[i]*thermo_.WiValue(speciei);
95  }
96 
97  forAll(reaction_.rhs(), i)
98  {
99  const label speciei = reaction_.rhs()[i].index;
100  Yprod0_[speciei] = thermo_.WiValue(speciei)/Wm*Xi[i];
101  }
102 
103  Info << "Maximum products mass concentrations: " << nl;
104  forAll(Yprod0_, i)
105  {
106  if (Yprod0_[i] > 0)
107  {
108  Info<< " " << thermo_.species()[i] << ": " << Yprod0_[i] << nl;
109  }
110  }
111 
112  // Normalise the stoichiometric coeff to mass
113  forAll(specieStoichCoeffs_, i)
114  {
115  specieStoichCoeffs_[i] =
116  specieStoichCoeffs_[i]*thermo_.WiValue(i)
117  /(thermo_.WiValue(fuelIndex_)*mag(specieStoichCoeffs_[fuelIndex_]));
118  }
119 }
120 
121 
122 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
123 
125 (
126  const word& modelType,
129  const word& combustionProperties
130 )
131 :
132  combustionModel(modelType, thermo, turb, combustionProperties),
133  reaction_(thermo_.species(), this->subDict("reaction")),
134  stoicRatio_(dimensionedScalar("stoicRatio", dimless, 0)),
135  s_(dimensionedScalar("s", dimless, 0)),
136  qFuel_(dimensionedScalar("qFuel", sqr(dimVelocity), 0)),
137  specieStoichCoeffs_(thermo_.species().size(), 0.0),
138  Yprod0_(thermo_.species().size(), 0.0),
139  fres_(Yprod0_.size()),
140  fuelIndex_(thermo_.species()[thermo.properties().lookup("fuel")]),
141  specieProd_(Yprod0_.size(), 1),
142  wFuel_
143  (
144  IOobject
145  (
146  this->thermo().phasePropertyName("wFuel"),
147  this->mesh().time().name(),
148  this->mesh(),
149  IOobject::NO_READ,
150  IOobject::NO_WRITE
151  ),
152  this->mesh(),
154  ),
155  semiImplicit_(readBool(this->coeffs_.lookup("semiImplicit")))
156 {
157  forAll(fres_, fresI)
158  {
159  IOobject header
160  (
161  "fres_" + thermo_.species()[fresI],
162  this->mesh().time().name(),
163  this->mesh()
164  );
165 
166  fres_.set
167  (
168  fresI,
169  new volScalarField
170  (
171  header,
172  this->mesh(),
173  dimensionedScalar("fres" + Foam::name(fresI), dimless, 0)
174  )
175  );
176  }
177 
178  calculateqFuel();
179 
181 
183 
184  if (semiImplicit_)
185  {
186  Info<< "Combustion mode: semi-implicit" << endl;
187  }
188  else
189  {
190  Info<< "Combustion mode: explicit" << endl;
191  }
192 }
193 
194 
195 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
196 
198 {}
199 
200 
201 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
202 
205 {
206  return wFuel_()*specieStoichCoeffs()[speciei];
207 }
208 
209 
212 {
213  const label specieI = thermo_.species()[Y.member()];
214 
215  volScalarField wSpecie
216  (
217  wFuel_*specieStoichCoeffs()[specieI]
218  );
219 
220  if (semiImplicit_)
221  {
222  const label fNorm = specieProd()[specieI];
223  const volScalarField fres(this->fres(specieI));
224  wSpecie /= max(fNorm*(Y - fres), scalar(1e-2));
225 
226  return -fNorm*wSpecie*fres + fNorm*fvm::Sp(wSpecie, Y);
227  }
228  else
229  {
230  return wSpecie + fvm::Sp(0.0*wSpecie, Y);
231  }
232 }
233 
234 
237 {
238  const label fuelI = fuelIndex();
239  volScalarField& YFuel =
240  const_cast<volScalarField&>(this->thermo().Y(fuelI));
241 
242  return -qFuel()*(R(YFuel) & YFuel);
243 }
244 
245 
247 {
248  if (combustionModel::read())
249  {
250  return true;
251  }
252  else
253  {
254  return false;
255  }
256 }
257 
258 
260 {
261  const label O2Index = thermo_.species()["O2"];
262  const volScalarField& YFuel = thermo_.Y()[fuelIndex_];
263  const volScalarField& YO2 = thermo_.Y()[O2Index];
264 
265  // reactants
266  forAll(reaction_.lhs(), i)
267  {
268  const label speciei = reaction_.lhs()[i].index;
269  if (speciei == fuelIndex_)
270  {
271  fres_[speciei] = max(YFuel - YO2/s_, scalar(0));
272  }
273  else if (speciei == O2Index)
274  {
275  fres_[speciei] = max(YO2 - YFuel*s_, scalar(0));
276  }
277  }
278 
279  // products
280  forAll(reaction_.rhs(), i)
281  {
282  const label speciei = reaction_.rhs()[i].index;
283  if (speciei != thermo_.defaultSpecie())
284  {
285  forAll(fres_[speciei], celli)
286  {
287  if (fres_[fuelIndex_][celli] > 0.0)
288  {
289  // rich mixture
290  fres_[speciei][celli] =
291  Yprod0_[speciei]
292  * (1.0 + YO2[celli]/s_.value() - YFuel[celli]);
293  }
294  else
295  {
296  // lean mixture
297  fres_[speciei][celli] =
298  Yprod0_[speciei]
299  * (
300  1.0
301  - YO2[celli]/s_.value()*stoicRatio_.value()
302  + YFuel[celli]*stoicRatio_.value()
303  );
304  }
305  }
306  }
307  }
308 }
309 
310 
311 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:433
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
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
const dimensionSet dimless
messageStream Info
void mag(LagrangianPatchField< scalar > &f, const LagrangianPatchField< Type > &f1)
const dimensionSet dimTime
const dimensionSet dimVolume
layerAndWeight max(const layerAndWeight &a, const layerAndWeight &b)
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
void sqr(LagrangianPatchField< typename outerProduct< Type, Type >::type > &f, const LagrangianPatchField< Type > &f1)
word name(const LagrangianState state)
Return a string representation of a Lagrangian state enumeration.
const dimensionSet dimMass
const dimensionSet dimVelocity
static const char nl
Definition: Ostream.H:267
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
PtrList< volScalarField > & Y
fluidMulticomponentThermo & thermo
Definition: createFields.H:31