diffusionMulticomponent.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) 2015-2018 OpenCFD Ltd.
9 -------------------------------------------------------------------------------
10 License
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 "fvcGrad.H"
30 #include "reactingMixture.H"
31 #include "fvCFD.H"
32 #include "mathematicalConstants.H"
33 
34 // * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
35 
36 template<class ReactionThermo, class ThermoType>
39 {
40  // Load default values
41  this->coeffs().readIfPresent("Ci", Ci_);
42  this->coeffs().readIfPresent("YoxStream", YoxStream_);
43  this->coeffs().readIfPresent("YfStream", YfStream_);
44  this->coeffs().readIfPresent("sigma", sigma_);
45  this->coeffs().readIfPresent("ftCorr", ftCorr_);
46  this->coeffs().readIfPresent("alpha", alpha_);
47  this->coeffs().readIfPresent("laminarIgn", laminarIgn_);
48 
49  typedef typename Reaction<ThermoType>::specieCoeffs specieCoeffs;
50 
51  const speciesTable& species = this->thermo().composition().species();
52 
53  scalarList specieStoichCoeffs(species.size());
54  const label nReactions = reactions_.size();
55 
56  for (label k=0; k < nReactions; k++)
57  {
58  RijPtr_.set
59  (
60  k,
61  new volScalarField
62  (
63  this->mesh_.newIOobject("Rijk" + Foam::name(k)),
64  this->mesh_,
67  )
68  );
69 
70  RijPtr_[k].storePrevIter();
71 
72  const List<specieCoeffs>& lhs = reactions_[k].lhs();
73  const List<specieCoeffs>& rhs = reactions_[k].rhs();
74 
75  const label fuelIndex = species.find(fuelNames_[k]);
76  const label oxidantIndex = species.find(oxidantNames_[k]);
77 
78  const scalar Wu = specieThermo_[fuelIndex].W();
79  const scalar Wox = specieThermo_[oxidantIndex].W();
80 
81  forAll(lhs, i)
82  {
83  const label specieI = lhs[i].index;
84  specieStoichCoeffs[specieI] = -lhs[i].stoichCoeff;
85  qFuel_[k] +=
86  specieThermo_[specieI].hc()*lhs[i].stoichCoeff/Wu;
87  }
88 
89  forAll(rhs, i)
90  {
91  const label specieI = rhs[i].index;
92  specieStoichCoeffs[specieI] = rhs[i].stoichCoeff;
93  qFuel_[k] -=
94  specieThermo_[specieI].hc()*rhs[i].stoichCoeff/Wu;
95  }
96 
97  Info << "Fuel heat of combustion : " << qFuel_[k] << endl;
98 
99  s_[k] =
100  (Wox*mag(specieStoichCoeffs[oxidantIndex]))
101  / (Wu*mag(specieStoichCoeffs[fuelIndex]));
102 
103  Info << "stoichiometric oxygen-fuel ratio : " << s_[k] << endl;
104 
105  stoicRatio_[k] = s_[k]*YfStream_[k]/YoxStream_[k];
106 
107  Info << "stoichiometric air-fuel ratio : " << stoicRatio_[k] << endl;
108 
109  const scalar fStoich = 1.0/(1.0 + stoicRatio_[k]);
110 
111  Info << "stoichiometric mixture fraction : " << fStoich << endl;
112  }
113 }
115 
116 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
117 
118 template<class ReactionThermo, class ThermoType>
121 (
122  const word& modelType,
123  ReactionThermo& thermo,
125  const word& combustionProperties
126 )
127 :
128  ChemistryCombustion<ReactionThermo>
129  (
130  modelType,
131  thermo,
132  turb,
133  combustionProperties
134  ),
135  reactions_
136  (
137  dynamic_cast<const reactingMixture<ThermoType>&>(thermo)
138  ),
139  specieThermo_
140  (
141  dynamic_cast<const reactingMixture<ThermoType>&>(thermo).
142  speciesData()
143  ),
144  RijPtr_(reactions_.size()),
145  Ci_(reactions_.size(), 1.0),
146  fuelNames_(this->coeffs().lookup("fuels")),
147  oxidantNames_(this->coeffs().lookup("oxidants")),
148  qFuel_(reactions_.size()),
149  stoicRatio_(reactions_.size()),
150  s_(reactions_.size()),
151  YoxStream_(reactions_.size(), 0.23),
152  YfStream_(reactions_.size(), 1.0),
153  sigma_(reactions_.size(), 0.02),
154  oxidantRes_(this->coeffs().lookup("oxidantRes")),
155  ftCorr_(reactions_.size(), Zero),
156  alpha_(1),
157  laminarIgn_(false)
158 {
159  init();
160 }
161 
162 
163 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
164 
165 template<class ReactionThermo, class ThermoType>
168 {
169  return this->chemistryPtr_->tc();
170 }
171 
172 
173 template<class ReactionThermo, class ThermoType>
176 {
177  if (this->active())
178  {
179  typedef typename Reaction<ThermoType>::specieCoeffs specieCoeffs;
180  const speciesTable& species = this->thermo().composition().species();
181 
182  const label nReactions = reactions_.size();
183 
184  PtrList<volScalarField> RijlPtr(nReactions);
185 
186  for (label k=0; k < nReactions; k++)
187  {
188  RijlPtr.set
189  (
190  k,
191  new volScalarField
192  (
193  this->mesh_.newIOobject("Rijl" + Foam::name(k)),
194  this->mesh_,
197  )
198  );
199 
200  volScalarField& Rijl = RijlPtr[k];
201 
202  // Obtain individual reaction rates for each reaction
203  const label fuelIndex = species.find(fuelNames_[k]);
204 
205  if (laminarIgn_)
206  {
207  Rijl.ref() = -this->chemistryPtr_->calculateRR(k, fuelIndex);
208  }
209 
210 
211  // Look for the fuelStoic
212  const List<specieCoeffs>& rhs = reactions_[k].rhs();
213  const List<specieCoeffs>& lhs = reactions_[k].lhs();
214 
215  // Set to zero RR's
216  forAll(lhs, l)
217  {
218  const label lIndex = lhs[l].index;
219  this->chemistryPtr_->RR(lIndex) =
221  }
222 
223  forAll(rhs, l)
224  {
225  const label rIndex = rhs[l].index;
226  this->chemistryPtr_->RR(rIndex) =
228  }
229  }
230 
231 
232  for (label k=0; k < nReactions; k++)
233  {
234  const label fuelIndex = species.find(fuelNames_[k]);
235  const label oxidantIndex = species.find(oxidantNames_[k]);
236 
237  const volScalarField& Yfuel =
238  this->thermo().composition().Y(fuelIndex);
239 
240  const volScalarField& Yox =
241  this->thermo().composition().Y(oxidantIndex);
242 
243  const volScalarField ft
244  (
245  "ft" + Foam::name(k),
246  (
247  s_[k]*Yfuel - (Yox - YoxStream_[k])
248  )
249  /(
250  s_[k]*YfStream_[k] + YoxStream_[k]
251  )
252  );
253 
254  const scalar sigma = sigma_[k];
255 
256  const scalar fStoich = 1.0/(1.0 + stoicRatio_[k]) + ftCorr_[k];
257 
258  const volScalarField OAvailScaled
259  (
260  "OAvailScaled",
261  Yox/max(oxidantRes_[k], 1e-3)
262  );
263 
264  const volScalarField preExp
265  (
266  "preExp" + Foam::name(k),
267  1.0 + sqr(OAvailScaled)
268  );
269 
270  const volScalarField filter
271  (
273  * exp(-sqr(ft - fStoich)/(2*sqr(sigma)))
274  );
275 
276  const volScalarField topHatFilter(pos(filter - 1e-3));
277 
278  const volScalarField prob("prob" + Foam::name(k), preExp*filter);
279 
280  const volScalarField RijkDiff
281  (
282  "RijkDiff",
283  Ci_[k]*this->turbulence().muEff()*prob*
284  (
285  mag(fvc::grad(Yfuel) & fvc::grad(Yox))
286  )
287  *pos(Yox)*pos(Yfuel)
288  );
289 
290  volScalarField& Rijk = RijPtr_[k];
291 
292  if (laminarIgn_)
293  {
294  Rijk =
295  min(RijkDiff, topHatFilter*RijlPtr[k]*pos(Yox)*pos(Yfuel));
296  }
297  else
298  {
299  Rijk = RijkDiff;
300  }
301 
302  Rijk.relax(alpha_);
303 
304  if (debug && this->mesh_.time().writeTime())
305  {
306  Rijk.write();
307  ft.write();
308  }
309 
310  // Look for the fuelStoic
311  const List<specieCoeffs>& rhs = reactions_[k].rhs();
312  const List<specieCoeffs>& lhs = reactions_[k].lhs();
313 
314  scalar fuelStoic = 1.0;
315  forAll(lhs, l)
316  {
317  const label lIndex = lhs[l].index;
318  if (lIndex == fuelIndex)
319  {
320  fuelStoic = lhs[l].stoichCoeff;
321  break;
322  }
323  }
324 
325  const scalar MwFuel = specieThermo_[fuelIndex].W();
326 
327  // Update left hand side species
328  forAll(lhs, l)
329  {
330  const label lIndex = lhs[l].index;
331 
332  const scalar stoichCoeff = lhs[l].stoichCoeff;
333 
334  this->chemistryPtr_->RR(lIndex) +=
335  -Rijk*stoichCoeff*specieThermo_[lIndex].W()/fuelStoic/MwFuel;
336 
337  }
338 
339  // Update right hand side species
340  forAll(rhs, r)
341  {
342  const label rIndex = rhs[r].index;
343 
344  const scalar stoichCoeff = rhs[r].stoichCoeff;
345 
346  this->chemistryPtr_->RR(rIndex) +=
347  Rijk*stoichCoeff*specieThermo_[rIndex].W()/fuelStoic/MwFuel;
348  }
349  }
350  }
351 }
352 
353 
354 template<class ReactionThermo, class ThermoType>
357 (
359 ) const
360 {
362  auto& Su = tSu.ref();
363 
364  if (this->active())
365  {
366  const label specieI =
367  this->thermo().composition().species().find(Y.member());
368 
369  Su += this->chemistryPtr_->RR(specieI);
370  }
371 
372  return tSu;
373 }
374 
375 
376 template<class ReactionThermo, class ThermoType>
380 {
381  auto tQdot = volScalarField::New
382  (
383  "Qdot",
385  this->mesh(),
388  );
389 
390  if (this->active())
391  {
392  volScalarField& Qdot = tQdot.ref();
393  Qdot = this->chemistryPtr_->Qdot();
394  }
396  return tQdot;
397 }
398 
399 
400 template<class ReactionThermo, class ThermoType>
403 {
405  {
406  this->coeffs().readIfPresent("Ci", Ci_);
407  this->coeffs().readIfPresent("sigma", sigma_);
408  this->coeffs().readIfPresent("oxidantRes", oxidantRes_);
409  this->coeffs().readIfPresent("ftCorr", ftCorr_);
410  this->coeffs().readIfPresent("alpha", alpha_);
411  this->coeffs().readIfPresent("laminarIgn", laminarIgn_);
412  return true;
413  }
414 
415  return false;
416 }
417 
418 
419 // ************************************************************************* //
static const word & zeroGradientType() noexcept
The type name for zeroGradient patch fields.
Definition: fvPatchField.H:221
virtual void correct()
Correct combustion rate.
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:47
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:116
Foam::reactingMixture.
const word zeroGradientType
A zeroGradient patch field type.
virtual tmp< volScalarField > Qdot() const
Heat release rate calculated from fuel consumption rate matrix.
zeroField Su
Definition: alphaSuSp.H:1
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:40
dimensionedSymmTensor sqr(const dimensionedVector &dv)
compressible::turbulenceModel & turb
dimensionedScalar sqrt(const dimensionedScalar &ds)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:531
const speciesTable & species() const
Return the table of species.
label k
Boltzmann constant.
Lookup type of boundary radiation properties.
Definition: lookup.H:57
PtrList< volScalarField > & Y()
Return the mass-fraction fields.
dimensionedScalar sigma("sigma", dimMass/sqr(dimTime), transportProperties)
Hold specie index and its coefficients in the reaction rate expression.
Definition: Reaction.H:86
scalar Qdot
Definition: solveChemistry.H:2
static word member(const word &name)
Return member (name without the extension)
Definition: IOobject.C:313
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:421
psiReactionThermo & thermo
Definition: createFields.H:28
dimensionedScalar pos(const dimensionedScalar &ds)
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
const dimensionedScalar e
Elementary charge.
Definition: createFields.H:11
dynamicFvMesh & mesh
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for expressions::valueTypeCode::INVALID.
Definition: exprTraits.C:127
dimensionedScalar exp(const dimensionedScalar &ds)
Calculate the gradient of the given field.
A class for handling words, derived from Foam::string.
Definition: word.H:63
virtual bool write(const bool writeOnProc=true) const
Write using setting from DB.
static tmp< T > New(Args &&... args)
Construct tmp with forwarding arguments.
Definition: tmp.H:206
constexpr scalar pi(M_PI)
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:26
Chemistry model wrapper for combustion models.
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
static tmp< GeometricField< scalar, fvPatchField, volMesh > > New(const word &name, IOobjectOption::registerOption regOpt, const Mesh &mesh, const dimensionSet &dims, const word &patchFieldType=fvPatchField< scalar >::calculatedType())
Return tmp field (NO_READ, NO_WRITE) from name, mesh, dimensions and patch type. [Takes current timeN...
Abstract base class for turbulence models (RAS, LES and laminar).
int debug
Static debugging option.
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
Find an entry if present, and assign to T val. FatalIOError if it is found and the number of tokens i...
Diffusion based turbulent combustion model for multicomponent species.
Info<< "Reading field U\"<< endl;volVectorField U(IOobject("U", runTime.timeName(), mesh, IOobject::MUST_READ, IOobject::AUTO_WRITE), mesh);volScalarField rho(IOobject("rho", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::AUTO_WRITE), thermo.rho());volVectorField rhoU(IOobject("rhoU", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *U);volScalarField rhoE(IOobject("rhoE", runTime.timeName(), mesh, IOobject::NO_READ, IOobject::NO_WRITE), rho *(e+0.5 *magSqr(U)));surfaceScalarField pos(IOobject("pos", runTime.timeName(), mesh), mesh, dimensionedScalar("pos", dimless, 1.0));surfaceScalarField neg(IOobject("neg", runTime.timeName(), mesh), mesh, dimensionedScalar("neg", dimless, -1.0));surfaceScalarField phi("phi", fvc::flux(rhoU));Info<< "Creating turbulence model\"<< endl;autoPtr< compressible::turbulenceModel > turbulence(compressible::turbulenceModel::New(rho, U, phi, thermo))
Definition: createFields.H:94
const dimensionSet dimEnergy
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
void relax(const scalar alpha)
Relax field (for steady-state solution).
A wordList with hashed named lookup, which can be faster in some situations than using the normal lis...
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
virtual bool read()
Update properties from given dictionary.
A list of pointers to objects of type <T>, with allocation/deallocation management of the pointers...
Definition: List.H:55
PtrList< volScalarField > & Y
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
messageStream Info
Information stream (stdout output on master, null elsewhere)
Internal & ref(const bool updateAccessTime=true)
Same as internalFieldRef()
virtual basicSpecieMixture & composition()=0
Return the composition of the multi-component mixture.
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
A class for managing temporary objects.
Definition: HashPtrTable.H:50
label find(const word &val) const
Find index of the value (searches the hash).
Do not request registration (bool: false)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:127