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 | Copyright (C) 2015 OpenCFD Ltd
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 "fvcGrad.H"
28 #include "reactingMixture.H"
29 #include "fvCFD.H"
30 #include "mathematicalConstants.H"
31 
32 // * * * * * * * * * * * * * * Private Member Functions * * * * * * * * * * //
33 
34 template<class CombThermoType, class ThermoType>
37 {
38  // Load default values
39  if (this->coeffs().found("Ci"))
40  {
41  this->coeffs().lookup("Ci") >> Ci_;
42  }
43 
44  if (this->coeffs().found("YoxStream"))
45  {
46  this->coeffs().lookup("YoxStream") >> YoxStream_;
47  }
48 
49  if (this->coeffs().found("YfStream"))
50  {
51  this->coeffs().lookup("YfStream") >> YfStream_;
52  }
53 
54  if (this->coeffs().found("sigma"))
55  {
56  this->coeffs().lookup("sigma") >> sigma_;
57  }
58 
59  if (this->coeffs().found("ftCorr"))
60  {
61  this->coeffs().lookup("ftCorr") >> ftCorr_;
62  }
63 
64  if (this->coeffs().found("alpha"))
65  {
66  alpha_ = readScalar(this->coeffs().lookup("alpha"));
67  }
68 
69  if (this->coeffs().found("laminarIgn"))
70  {
71  this->coeffs().lookup("laminarIgn") >> laminarIgn_;
72  }
73 
74  typedef typename Reaction<ThermoType>::specieCoeffs specieCoeffs;
75 
76  const speciesTable& species = this->thermo().composition().species();
77 
78  scalarList specieStoichCoeffs(species.size());
79  const label nReactions = reactions_.size();
80 
81  for (label k=0; k < nReactions; k++)
82  {
83  RijPtr_.set
84  (
85  k,
86  new volScalarField
87  (
88  IOobject
89  (
90  "Rijk" + name(k),
91  this->mesh_.time().timeName(),
92  this->mesh_,
93  IOobject::NO_READ,
94  IOobject::NO_WRITE,
95  false
96  ),
97  this->mesh_,
99  zeroGradientFvPatchScalarField::typeName
100  )
101  );
102 
103  RijPtr_[k].storePrevIter();
104 
105  const List<specieCoeffs>& lhs = reactions_[k].lhs();
106  const List<specieCoeffs>& rhs = reactions_[k].rhs();
107 
108  const label fuelIndex = species[fuelNames_[k]];
109  const label oxidantIndex = species[oxidantNames_[k]];
110 
111  const scalar Wu = specieThermo_[fuelIndex].W();
112  const scalar Wox = specieThermo_[oxidantIndex].W();
113 
114  forAll(lhs, i)
115  {
116  const label specieI = lhs[i].index;
117  specieStoichCoeffs[specieI] = -lhs[i].stoichCoeff;
118  qFuel_[k] +=
119  specieThermo_[specieI].hc()*lhs[i].stoichCoeff/Wu;
120  }
121 
122  forAll(rhs, i)
123  {
124  const label specieI = rhs[i].index;
125  specieStoichCoeffs[specieI] = rhs[i].stoichCoeff;
126  qFuel_[k] -=
127  specieThermo_[specieI].hc()*rhs[i].stoichCoeff/Wu;
128  }
129 
130  Info << "Fuel heat of combustion : " << qFuel_[k] << endl;
131 
132  s_[k] =
133  (Wox*mag(specieStoichCoeffs[oxidantIndex]))
134  / (Wu*mag(specieStoichCoeffs[fuelIndex]));
135 
136  Info << "stoichiometric oxygen-fuel ratio : " << s_[k] << endl;
137 
138  stoicRatio_[k] = s_[k]*YfStream_[k]/YoxStream_[k];
139 
140  Info << "stoichiometric air-fuel ratio : " << stoicRatio_[k] << endl;
141 
142  const scalar fStoich = 1.0/(1.0 + stoicRatio_[k]);
143 
144  Info << "stoichiometric mixture fraction : " << fStoich << endl;
145 
146  }
147 }
148 
149 
150 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
151 
152 template<class CombThermoType, class ThermoType>
155 (
156  const word& modelType, const fvMesh& mesh, const word& phaseName
157 )
158 :
159  CombThermoType(modelType, mesh, phaseName),
160  reactions_
161  (
162  dynamic_cast<const reactingMixture<ThermoType>&>(this->thermo())
163  ),
164  specieThermo_
165  (
166  dynamic_cast<const reactingMixture<ThermoType>&>
167  (this->thermo()).speciesData()
168  ),
169  RijPtr_(reactions_.size()),
170  Ci_(reactions_.size(), 1.0),
171  fuelNames_(this->coeffs().lookup("fuels")),
172  oxidantNames_(this->coeffs().lookup("oxidants")),
173  qFuel_(reactions_.size()),
174  stoicRatio_(reactions_.size()),
175  s_(reactions_.size()),
176  YoxStream_(reactions_.size(), 0.23),
177  YfStream_(reactions_.size(), 1.0),
178  sigma_(reactions_.size(), 0.02),
179  oxidantRes_(this->coeffs().lookup("oxidantRes")),
180  ftCorr_(reactions_.size(), 0.0),
181  alpha_(1),
182  laminarIgn_(false)
183 {
184  init();
185 }
186 
187 
188 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
189 
190 template<class CombThermoType, class ThermoType>
193 {}
194 
195 
196 // * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * * //
197 
198 template<class CombThermoType, class ThermoType>
201 {
202  return this->chemistryPtr_->tc();
203 }
204 
205 
206 template<class CombThermoType, class ThermoType>
209 {
210  if (this->active())
211  {
212  typedef typename Reaction<ThermoType>::specieCoeffs specieCoeffs;
213  const speciesTable& species = this->thermo().composition().species();
214 
215  const label nReactions = reactions_.size();
216 
217  PtrList<volScalarField> RijlPtr(nReactions);
218 
219  for (label k=0; k < nReactions; k++)
220  {
221 
222  RijlPtr.set
223  (
224  k,
225  new volScalarField
226  (
227  IOobject
228  (
229  "Rijl" + word(k),
230  this->mesh_.time().timeName(),
231  this->mesh_,
232  IOobject::NO_READ,
233  IOobject::NO_WRITE,
234  false
235  ),
236  this->mesh_,
238  zeroGradientFvPatchScalarField::typeName
239  )
240  );
241 
242  volScalarField& Rijl = RijlPtr[k];
243 
244  // Obtain individual reaction rates for each reaction
245  const label fuelIndex = species[fuelNames_[k]];
246 
247  if (laminarIgn_)
248  {
249  Rijl.dimensionedInternalField() =
250  -this->chemistryPtr_->calculateRR(k, fuelIndex);
251  }
252 
253 
254  // Look for the fuelStoic
255  const List<specieCoeffs>& rhs = reactions_[k].rhs();
256  const List<specieCoeffs>& lhs = reactions_[k].lhs();
257 
258  // Set to zero RR's
259  forAll (lhs, l)
260  {
261  const label lIndex = lhs[l].index;
262  this->chemistryPtr_->RR(lIndex) =
264  }
265 
266  forAll (rhs, l)
267  {
268  const label rIndex = rhs[l].index;
269  this->chemistryPtr_->RR(rIndex) =
271  }
272  }
273 
274 
275  for (label k=0; k < nReactions; k++)
276  {
277  const label fuelIndex = species[fuelNames_[k]];
278  const label oxidantIndex = species[oxidantNames_[k]];
279 
280  const volScalarField& Yfuel =
281  this->thermo().composition().Y(fuelIndex);
282 
283  const volScalarField& Yox =
284  this->thermo().composition().Y(oxidantIndex);
285 
286  const volScalarField ft
287  (
288  "ft" + name(k),
289  (
290  s_[k]*Yfuel
291  - (Yox - YoxStream_[k])
292  )
293  /
294  (
295  s_[k]*YfStream_[k] + YoxStream_[k]
296  )
297  );
298 
299  const scalar sigma = sigma_[k];
300 
301  const scalar fStoich = 1.0/(1.0 + stoicRatio_[k]) + ftCorr_[k];
302 
303  const volScalarField OAvailScaled
304  (
305  "OAvailScaled",
306  Yox/max(oxidantRes_[k], 1e-3)
307  );
308 
309  const volScalarField preExp
310  (
311  "preExp" + name(k),
312  1.0 + sqr(OAvailScaled)
313  );
314 
315  const volScalarField filter
316  (
318  * exp(-sqr(ft - fStoich)/(2*sqr(sigma)))
319  );
320 
321  const volScalarField topHatFilter
322  (
323  pos(filter - 1e-3)
324  );
325 
326  const volScalarField prob
327  (
328  "prob" + name(k), preExp*filter
329  );
330 
331  const volScalarField RijkDiff
332  (
333  "RijkDiff",
334  Ci_[k]*this->turbulence().muEff()*prob*
335  (
336  mag(fvc::grad(Yfuel) & fvc::grad(Yox))
337  )
338  *pos(Yox)*pos(Yfuel)
339  );
340 
341  volScalarField& Rijk = RijPtr_[k];
342 
343  if (laminarIgn_)
344  {
345  Rijk =
346  min(RijkDiff, topHatFilter*RijlPtr[k]*pos(Yox)*pos(Yfuel));
347  }
348  else
349  {
350  Rijk = RijkDiff;
351  }
352 
353  Rijk.relax(alpha_);
354 
355  if (this->mesh_.time().outputTime() && debug)
356  {
357  Rijk.write();
358  ft.write();
359  }
360 
361  // Look for the fuelStoic
362  const List<specieCoeffs>& rhs = reactions_[k].rhs();
363  const List<specieCoeffs>& lhs = reactions_[k].lhs();
364 
365  scalar fuelStoic = 1.0;
366  forAll (lhs, l)
367  {
368  const label lIndex = lhs[l].index;
369  if (lIndex == fuelIndex)
370  {
371  fuelStoic = lhs[l].stoichCoeff;
372  break;
373  }
374  }
375 
376  const scalar MwFuel = specieThermo_[fuelIndex].W();
377 
378  // Update left hand side species
379  forAll (lhs, l)
380  {
381  const label lIndex = lhs[l].index;
382 
383  const scalar stoichCoeff = lhs[l].stoichCoeff;
384 
385  this->chemistryPtr_->RR(lIndex) +=
386  -Rijk*stoichCoeff*specieThermo_[lIndex].W()
387  /fuelStoic/MwFuel;
388 
389  }
390 
391  // Update right hand side species
392  forAll (rhs, r)
393  {
394  const label rIndex = rhs[r].index;
395 
396  const scalar stoichCoeff = rhs[r].stoichCoeff;
397 
398  this->chemistryPtr_->RR(rIndex) +=
399  Rijk*stoichCoeff*specieThermo_[rIndex].W()
400  /fuelStoic/MwFuel;
401  }
402  }
403  }
404 }
405 
406 
407 template<class CombThermoType, class ThermoType>
411 {
413 
414  fvScalarMatrix& Su = tSu();
415 
416  if (this->active())
417  {
418  const label specieI = this->thermo().composition().species()[Y.name()];
419  Su += this->chemistryPtr_->RR(specieI);
420  }
421 
422  return tSu;
423 }
424 
425 
426 template<class CombThermoType, class ThermoType>
429 dQ() const
430 {
432  (
433  new volScalarField
434  (
435  IOobject
436  (
437  "dQ",
438  this->mesh().time().timeName(),
439  this->mesh(),
440  IOobject::NO_READ,
441  IOobject::NO_WRITE,
442  false
443  ),
444  this->mesh(),
445  dimensionedScalar("dQ", dimEnergy/dimTime, 0.0),
446  zeroGradientFvPatchScalarField::typeName
447  )
448  );
449 
450  if (this->active())
451  {
452  volScalarField& dQ = tdQ();
453  dQ = this->chemistryPtr_->dQ();
454  }
455 
456  return tdQ;
457 }
458 
459 
460 template<class CombThermoType, class ThermoType>
463 Sh() const
464 {
466  (
467  new volScalarField
468  (
469  IOobject
470  (
471  "Sh",
472  this->mesh().time().timeName(),
473  this->mesh(),
474  IOobject::NO_READ,
475  IOobject::NO_WRITE,
476  false
477  ),
478  this->mesh(),
480  zeroGradientFvPatchScalarField::typeName
481  )
482  );
483 
484  if (this->active())
485  {
486  scalarField& Sh = tSh();
487  Sh = this->chemistryPtr_->Sh();
488  }
489 
490  return tSh;
491 }
492 
493 
494 template<class CombThermoType, class ThermoType>
497 {
498  if (CombThermoType::read())
499  {
500  this->coeffs().readIfPresent("Ci", Ci_);
501  this->coeffs().readIfPresent("sigma", sigma_);
502  this->coeffs().readIfPresent("oxidantRes", oxidantRes_);
503  this->coeffs().readIfPresent("ftCorr", ftCorr_);
504  this->coeffs().readIfPresent("alpha", alpha_);
505  this->coeffs().readIfPresent("laminarIgn", laminarIgn_);
506  return true;
507  }
508  else
509  {
510  return false;
511  }
512 }
513 
514 
515 // ************************************************************************* //
Foam::IOobject
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:91
mathematicalConstants.H
Foam::combustionModels::diffusionMulticomponent::init
void init()
Initialize.
Definition: diffusionMulticomponent.C:36
Foam::combustionModels::diffusionMulticomponent::~diffusionMulticomponent
virtual ~diffusionMulticomponent()
Destructor.
Definition: diffusionMulticomponent.C:192
Foam::fvc::grad
tmp< GeometricField< typename outerProduct< vector, Type >::type, fvPatchField, volMesh >> grad(const GeometricField< Type, fvsPatchField, surfaceMesh > &ssf)
Definition: fvcGrad.C:52
reactingMixture.H
Foam::GeometricField::dimensionedInternalField
DimensionedInternalField & dimensionedInternalField()
Return dimensioned internal field.
Definition: GeometricField.C:713
Foam::word
A class for handling words, derived from string.
Definition: word.H:59
Sh
scalar Sh
Definition: solveChemistry.H:2
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:406
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:118
Foam::dimEnergy
const dimensionSet dimEnergy
Foam::read
bool read(const char *, int32_t &)
Definition: int32IO.C:87
Foam::reactingMixture
Foam::reactingMixture.
Definition: reactingMixture.H:51
dQ
dQ
Definition: YEEqn.H:14
Foam::Reaction::specieCoeffs
Class to hold the specie index and its coefficients in the.
Definition: Reaction.H:91
Foam::combustionModels::diffusionMulticomponent::R
virtual tmp< fvScalarMatrix > R(volScalarField &Y) const
Fuel consumption rate matrix.
Definition: diffusionMulticomponent.C:410
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:252
Foam::combustionModels::diffusionMulticomponent::diffusionMulticomponent
diffusionMulticomponent(const diffusionMulticomponent &)
Disallow copy construct.
Foam::mag
dimensioned< scalar > mag(const dimensioned< Type > &)
thermo
rhoThermo & thermo
Definition: setRegionFluidFields.H:3
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:51
Foam::PtrList::set
bool set(const label) const
Is element set.
diffusionMulticomponent.H
Foam::label
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
Foam::fvScalarMatrix
fvMatrix< scalar > fvScalarMatrix
Definition: fvMatricesFwd.H:42
Foam::Field
Pre-declare SubField and related Field type.
Definition: Field.H:57
Foam::fvc::Su
tmp< GeometricField< Type, fvPatchField, volMesh > > Su(const GeometricField< Type, fvPatchField, volMesh > &su, const GeometricField< Type, fvPatchField, volMesh > &vf)
Definition: fvcSup.C:44
Foam::hashedWordList
A wordList with hashed indices for faster lookup by name.
Definition: hashedWordList.H:57
turbulence
autoPtr< compressible::turbulenceModel > turbulence
Definition: createFields.H:23
Foam::Info
messageStream Info
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:41
Foam::PtrList
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: List.H:61
Foam::combustionModels::diffusionMulticomponent::correct
virtual void correct()
Correct combustion rate.
Definition: diffusionMulticomponent.C:208
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:18
Foam::dimensioned
Generic dimensioned Type class.
Definition: dimensionedScalarFwd.H:41
Foam::fvMesh
Mesh data needed to do the Finite Volume discretisation.
Definition: fvMesh.H:78
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:49
Foam::combustionModels::diffusionMulticomponent::dQ
virtual tmp< volScalarField > dQ() const
Heat release rate calculated from fuel consumption rate matrix.
Definition: diffusionMulticomponent.C:429
Foam::combustionModels::diffusionMulticomponent::read
virtual bool read()
Update properties from given dictionary.
Definition: diffusionMulticomponent.C:496
readScalar
#define readScalar
Definition: doubleScalar.C:38
Foam::max
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::combustionModels::diffusionMulticomponent::Sh
virtual tmp< volScalarField > Sh() const
Return source for enthalpy equation [kg/m/s3].
Definition: diffusionMulticomponent.C:463
Foam::sqr
dimensionedSymmTensor sqr(const dimensionedVector &dv)
Definition: dimensionedSymmTensor.C:49
Foam::combustionModels::diffusionMulticomponent::tc
tmp< volScalarField > tc() const
Return the chemical time scale.
Definition: diffusionMulticomponent.C:200
Foam::List
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
Foam::sqrt
dimensionedScalar sqrt(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:142
fvcGrad.H
Calculate the gradient of the given field.
k
label k
Boltzmann constant.
Definition: LISASMDCalcMethod2.H:41
Foam::constant::electromagnetic::e
const dimensionedScalar e
Elementary charge.
Definition: doubleFloat.H:94
timeName
word timeName
Definition: getTimeIndex.H:3
Foam::GeometricField::relax
void relax(const scalar alpha)
Relax field (for steady-state solution).
Definition: GeometricField.C:917
Foam::fvMatrix
A special matrix type and solver, designed for finite volume solutions of scalar equations....
Definition: fvPatchField.H:68
Foam::constant::physicoChemical::sigma
const dimensionedScalar sigma
Stefan-Boltzmann constant: default SI units: [W/m2/K4].
Foam::constant::mathematical::pi
const scalar pi(M_PI)
fvCFD.H
Foam::List::size
void size(const label)
Override size to be inconsistent with allocated storage.
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:58
Foam::GeometricField
Generic GeometricField class.
Definition: surfaceFieldsFwd.H:52
Foam::min
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
Foam::name
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
Y
PtrList< volScalarField > & Y
Definition: createFields.H:36
lookup
stressControl lookup("compactNormalStress") >> compactNormalStress
Foam::pos
dimensionedScalar pos(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:190