StandardChemistryModel.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) 2011-2017 OpenFOAM Foundation
9  Copyright (C) 2020-2021 OpenCFD Ltd.
10 -------------------------------------------------------------------------------
11 License
12  This file is part of OpenFOAM.
13 
14  OpenFOAM is free software: you can redistribute it and/or modify it
15  under the terms of the GNU General Public License as published by
16  the Free Software Foundation, either version 3 of the License, or
17  (at your option) any later version.
18 
19  OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
20  ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
21  FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
22  for more details.
23 
24  You should have received a copy of the GNU General Public License
25  along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
26 
27 \*---------------------------------------------------------------------------*/
28 
29 #include "StandardChemistryModel.H"
30 #include "reactingMixture.H"
31 #include "UniformField.H"
33 
34 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
35 
36 template<class ReactionThermo, class ThermoType>
38 (
39  ReactionThermo& thermo
40 )
41 :
43  ODESystem(),
44  Y_(this->thermo().composition().Y()),
45  reactions_
46  (
47  dynamic_cast<const reactingMixture<ThermoType>&>(this->thermo())
48  ),
49  specieThermo_
50  (
51  dynamic_cast<const reactingMixture<ThermoType>&>
52  (this->thermo()).speciesData()
53  ),
54 
55  nSpecie_(Y_.size()),
56  nReaction_(reactions_.size()),
57  Treact_
58  (
60  (
61  "Treact",
62  0.0
63  )
64  ),
65  RR_(nSpecie_),
66  c_(nSpecie_),
67  dcdt_(nSpecie_)
68 {
69  // Create the fields for the chemistry sources
70  forAll(RR_, fieldi)
71  {
72  RR_.set
73  (
74  fieldi,
76  (
77  IOobject
78  (
79  "RR." + Y_[fieldi].name(),
80  this->mesh().time().timeName(),
81  this->mesh(),
82  IOobject::NO_READ,
83  IOobject::NO_WRITE
84  ),
85  this->mesh(),
87  )
88  );
89  }
90 
91  Info<< "StandardChemistryModel: Number of species = " << nSpecie_
92  << " and reactions = " << nReaction_ << endl;
93 }
94 
95 
96 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
97 
98 template<class ReactionThermo, class ThermoType>
101 {}
102 
103 
104 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
105 
106 template<class ReactionThermo, class ThermoType>
108 (
109  const scalarField& c,
110  const scalar T,
111  const scalar p,
112  scalarField& dcdt
113 ) const
114 {
115  scalar pf, cf, pr, cr;
116  label lRef, rRef;
117 
118  dcdt = Zero;
119 
120  forAll(reactions_, i)
121  {
122  const Reaction<ThermoType>& R = reactions_[i];
123 
124  scalar omegai = omega
125  (
126  R, c, T, p, pf, cf, lRef, pr, cr, rRef
127  );
128 
129  forAll(R.lhs(), s)
130  {
131  const label si = R.lhs()[s].index;
132  const scalar sl = R.lhs()[s].stoichCoeff;
133  dcdt[si] -= sl*omegai;
134  }
135 
136  forAll(R.rhs(), s)
137  {
138  const label si = R.rhs()[s].index;
139  const scalar sr = R.rhs()[s].stoichCoeff;
140  dcdt[si] += sr*omegai;
141  }
142  }
143 }
144 
145 
146 template<class ReactionThermo, class ThermoType>
148 (
149  const label index,
150  const scalarField& c,
151  const scalar T,
152  const scalar p,
153  scalar& pf,
154  scalar& cf,
155  label& lRef,
156  scalar& pr,
157  scalar& cr,
158  label& rRef
159 ) const
160 {
161  const Reaction<ThermoType>& R = reactions_[index];
162  scalar w = omega(R, c, T, p, pf, cf, lRef, pr, cr, rRef);
163  return(w);
164 }
165 
166 
167 template<class ReactionThermo, class ThermoType>
169 (
170  const Reaction<ThermoType>& R,
171  const scalarField& c,
172  const scalar T,
173  const scalar p,
174  scalar& pf,
175  scalar& cf,
176  label& lRef,
177  scalar& pr,
178  scalar& cr,
179  label& rRef
180 ) const
181 {
182  const scalar kf = R.kf(p, T, c);
183  const scalar kr = R.kr(kf, p, T, c);
184 
185  pf = 1.0;
186  pr = 1.0;
187 
188  const label Nl = R.lhs().size();
189  const label Nr = R.rhs().size();
190 
191  label slRef = 0;
192  lRef = R.lhs()[slRef].index;
193 
194  pf = kf;
195  for (label s = 1; s < Nl; s++)
196  {
197  const label si = R.lhs()[s].index;
198 
199  if (c[si] < c[lRef])
200  {
201  const scalar exp = R.lhs()[slRef].exponent;
202  pf *= pow(max(c[lRef], 0.0), exp);
203  lRef = si;
204  slRef = s;
205  }
206  else
207  {
208  const scalar exp = R.lhs()[s].exponent;
209  pf *= pow(max(c[si], 0.0), exp);
210  }
211  }
212  cf = max(c[lRef], 0.0);
213 
214  {
215  const scalar exp = R.lhs()[slRef].exponent;
216  if (exp < 1.0)
217  {
218  if (cf > SMALL)
219  {
220  pf *= pow(cf, exp - 1.0);
221  }
222  else
223  {
224  pf = 0.0;
225  }
226  }
227  else
228  {
229  pf *= pow(cf, exp - 1.0);
230  }
231  }
232 
233  label srRef = 0;
234  rRef = R.rhs()[srRef].index;
235 
236  // Find the matrix element and element position for the rhs
237  pr = kr;
238  for (label s = 1; s < Nr; s++)
239  {
240  const label si = R.rhs()[s].index;
241  if (c[si] < c[rRef])
242  {
243  const scalar exp = R.rhs()[srRef].exponent;
244  pr *= pow(max(c[rRef], 0.0), exp);
245  rRef = si;
246  srRef = s;
247  }
248  else
249  {
250  const scalar exp = R.rhs()[s].exponent;
251  pr *= pow(max(c[si], 0.0), exp);
252  }
253  }
254  cr = max(c[rRef], 0.0);
255 
256  {
257  const scalar exp = R.rhs()[srRef].exponent;
258  if (exp < 1.0)
259  {
260  if (cr>SMALL)
261  {
262  pr *= pow(cr, exp - 1.0);
263  }
264  else
265  {
266  pr = 0.0;
267  }
268  }
269  else
270  {
271  pr *= pow(cr, exp - 1.0);
272  }
273  }
274 
275  return pf*cf - pr*cr;
276 }
277 
278 
279 template<class ReactionThermo, class ThermoType>
281 (
282  const scalar time,
283  const scalarField& c,
284  scalarField& dcdt
285 ) const
286 {
287  const scalar T = c[nSpecie_];
288  const scalar p = c[nSpecie_ + 1];
289 
290  forAll(c_, i)
291  {
292  c_[i] = max(c[i], 0.0);
293  }
294 
295  omega(c_, T, p, dcdt);
296 
297  // Constant pressure
298  // dT/dt = ...
299  scalar rho = 0.0;
300  for (label i = 0; i < nSpecie_; i++)
301  {
302  const scalar W = specieThermo_[i].W();
303  rho += W*c_[i];
304  }
305  scalar cp = 0.0;
306  for (label i=0; i<nSpecie_; i++)
307  {
308  cp += c_[i]*specieThermo_[i].cp(p, T);
309  }
310  cp /= rho;
311 
312  scalar dT = 0.0;
313  for (label i = 0; i < nSpecie_; i++)
314  {
315  const scalar hi = specieThermo_[i].ha(p, T);
316  dT += hi*dcdt[i];
317  }
318  dT /= rho*cp;
319 
320  dcdt[nSpecie_] = -dT;
321 
322  // dp/dt = ...
323  dcdt[nSpecie_ + 1] = 0.0;
324 }
325 
326 
327 template<class ReactionThermo, class ThermoType>
329 (
330  const scalar t,
331  const scalarField& c,
332  scalarField& dcdt,
333  scalarSquareMatrix& dfdc
334 ) const
335 {
336  const scalar T = c[nSpecie_];
337  const scalar p = c[nSpecie_ + 1];
338 
339  forAll(c_, i)
340  {
341  c_[i] = max(c[i], 0.0);
342  }
343 
344  dfdc = Zero;
345 
346  // Length of the first argument must be nSpecie_
347  omega(c_, T, p, dcdt);
348 
349  forAll(reactions_, ri)
350  {
351  const Reaction<ThermoType>& R = reactions_[ri];
352 
353  const scalar kf0 = R.kf(p, T, c_);
354  const scalar kr0 = R.kr(kf0, p, T, c_);
355 
356  forAll(R.lhs(), j)
357  {
358  const label sj = R.lhs()[j].index;
359  scalar kf = kf0;
360  forAll(R.lhs(), i)
361  {
362  const label si = R.lhs()[i].index;
363  const scalar el = R.lhs()[i].exponent;
364  if (i == j)
365  {
366  if (el < 1.0)
367  {
368  if (c_[si] > SMALL)
369  {
370  kf *= el*pow(c_[si], el - 1.0);
371  }
372  else
373  {
374  kf = 0.0;
375  }
376  }
377  else
378  {
379  kf *= el*pow(c_[si], el - 1.0);
380  }
381  }
382  else
383  {
384  kf *= pow(c_[si], el);
385  }
386  }
387 
388  forAll(R.lhs(), i)
389  {
390  const label si = R.lhs()[i].index;
391  const scalar sl = R.lhs()[i].stoichCoeff;
392  dfdc(si, sj) -= sl*kf;
393  }
394  forAll(R.rhs(), i)
395  {
396  const label si = R.rhs()[i].index;
397  const scalar sr = R.rhs()[i].stoichCoeff;
398  dfdc(si, sj) += sr*kf;
399  }
400  }
401 
402  forAll(R.rhs(), j)
403  {
404  const label sj = R.rhs()[j].index;
405  scalar kr = kr0;
406  forAll(R.rhs(), i)
407  {
408  const label si = R.rhs()[i].index;
409  const scalar er = R.rhs()[i].exponent;
410  if (i == j)
411  {
412  if (er < 1.0)
413  {
414  if (c_[si] > SMALL)
415  {
416  kr *= er*pow(c_[si], er - 1.0);
417  }
418  else
419  {
420  kr = 0.0;
421  }
422  }
423  else
424  {
425  kr *= er*pow(c_[si], er - 1.0);
426  }
427  }
428  else
429  {
430  kr *= pow(c_[si], er);
431  }
432  }
433 
434  forAll(R.lhs(), i)
435  {
436  const label si = R.lhs()[i].index;
437  const scalar sl = R.lhs()[i].stoichCoeff;
438  dfdc(si, sj) += sl*kr;
439  }
440  forAll(R.rhs(), i)
441  {
442  const label si = R.rhs()[i].index;
443  const scalar sr = R.rhs()[i].stoichCoeff;
444  dfdc(si, sj) -= sr*kr;
445  }
446  }
447  }
448 
449  // Calculate the dcdT elements numerically
450  const scalar delta = 1.0e-3;
451 
452  omega(c_, T + delta, p, dcdt_);
453  for (label i=0; i<nSpecie_; i++)
454  {
455  dfdc(i, nSpecie_) = dcdt_[i];
456  }
457 
458  omega(c_, T - delta, p, dcdt_);
459  for (label i=0; i<nSpecie_; i++)
460  {
461  dfdc(i, nSpecie_) = 0.5*(dfdc(i, nSpecie_) - dcdt_[i])/delta;
462  }
463 
464  dfdc(nSpecie_, nSpecie_) = 0;
465  dfdc(nSpecie_ + 1, nSpecie_) = 0;
466 }
467 
468 
469 template<class ReactionThermo, class ThermoType>
472 {
474  (
475  new volScalarField
476  (
477  IOobject
478  (
479  "tc",
480  this->time().timeName(),
481  this->mesh(),
482  IOobject::NO_READ,
483  IOobject::NO_WRITE,
484  false
485  ),
486  this->mesh(),
487  dimensionedScalar("small", dimTime, SMALL),
488  extrapolatedCalculatedFvPatchScalarField::typeName
489  )
490  );
491 
492  scalarField& tc = ttc.ref();
493 
494  tmp<volScalarField> trho(this->thermo().rho());
495  const scalarField& rho = trho();
496 
497  const scalarField& T = this->thermo().T();
498  const scalarField& p = this->thermo().p();
499 
500  const label nReaction = reactions_.size();
501 
502  scalar pf, cf, pr, cr;
503  label lRef, rRef;
504 
505  if (this->chemistry_)
506  {
507  forAll(rho, celli)
508  {
509  const scalar rhoi = rho[celli];
510  const scalar Ti = T[celli];
511  const scalar pi = p[celli];
512 
513  scalar cSum = 0.0;
514 
515  for (label i=0; i<nSpecie_; i++)
516  {
517  c_[i] = rhoi*Y_[i][celli]/specieThermo_[i].W();
518  cSum += c_[i];
519  }
520 
521  forAll(reactions_, i)
522  {
523  const Reaction<ThermoType>& R = reactions_[i];
524 
525  omega(R, c_, Ti, pi, pf, cf, lRef, pr, cr, rRef);
526 
527  forAll(R.rhs(), s)
528  {
529  tc[celli] += R.rhs()[s].stoichCoeff*pf*cf;
530  }
531  }
532 
533  tc[celli] = nReaction*cSum/tc[celli];
534  }
535  }
536 
537  ttc.ref().correctBoundaryConditions();
538 
539  return ttc;
540 }
541 
542 
543 template<class ReactionThermo, class ThermoType>
546 {
547  tmp<volScalarField> tQdot
548  (
549  new volScalarField
550  (
551  IOobject
552  (
553  "Qdot",
554  this->mesh_.time().timeName(),
555  this->mesh_,
556  IOobject::NO_READ,
557  IOobject::NO_WRITE,
558  false
559  ),
560  this->mesh_,
562  )
563  );
564 
565  if (this->chemistry_)
566  {
567  scalarField& Qdot = tQdot.ref();
568 
569  forAll(Y_, i)
570  {
571  forAll(Qdot, celli)
572  {
573  const scalar hi = specieThermo_[i].Hc();
574  Qdot[celli] -= hi*RR_[i][celli];
575  }
576  }
577  }
578 
579  return tQdot;
580 }
581 
582 
583 template<class ReactionThermo, class ThermoType>
586 (
587  const label ri,
588  const label si
589 ) const
590 {
591  scalar pf, cf, pr, cr;
592  label lRef, rRef;
593 
595  (
597  (
598  IOobject
599  (
600  "RR",
601  this->mesh().time().timeName(),
602  this->mesh(),
603  IOobject::NO_READ,
604  IOobject::NO_WRITE
605  ),
606  this->mesh(),
608  )
609  );
610 
612 
613  tmp<volScalarField> trho(this->thermo().rho());
614  const scalarField& rho = trho();
615 
616  const scalarField& T = this->thermo().T();
617  const scalarField& p = this->thermo().p();
618 
619  forAll(rho, celli)
620  {
621  const scalar rhoi = rho[celli];
622  const scalar Ti = T[celli];
623  const scalar pi = p[celli];
624 
625  for (label i=0; i<nSpecie_; i++)
626  {
627  const scalar Yi = Y_[i][celli];
628  c_[i] = rhoi*Yi/specieThermo_[i].W();
629  }
630 
631  const scalar w = omegaI
632  (
633  ri,
634  c_,
635  Ti,
636  pi,
637  pf,
638  cf,
639  lRef,
640  pr,
641  cr,
642  rRef
643  );
644 
645  RR[celli] = w*specieThermo_[si].W();
646  }
647 
648  return tRR;
649 }
650 
651 
652 template<class ReactionThermo, class ThermoType>
654 {
655  if (!this->chemistry_)
656  {
657  return;
658  }
659 
660  tmp<volScalarField> trho(this->thermo().rho());
661  const scalarField& rho = trho();
662 
663  const scalarField& T = this->thermo().T();
664  const scalarField& p = this->thermo().p();
665 
666  forAll(rho, celli)
667  {
668  const scalar rhoi = rho[celli];
669  const scalar Ti = T[celli];
670  const scalar pi = p[celli];
671 
672  for (label i=0; i<nSpecie_; i++)
673  {
674  const scalar Yi = Y_[i][celli];
675  c_[i] = rhoi*Yi/specieThermo_[i].W();
676  }
677 
678  omega(c_, Ti, pi, dcdt_);
679 
680  for (label i=0; i<nSpecie_; i++)
681  {
682  RR_[i][celli] = dcdt_[i]*specieThermo_[i].W();
683  }
684  }
685 }
686 
687 
688 template<class ReactionThermo, class ThermoType>
689 template<class DeltaTType>
691 (
692  const DeltaTType& deltaT
693 )
694 {
696 
697  scalar deltaTMin = GREAT;
698 
699  if (!this->chemistry_)
700  {
701  return deltaTMin;
702  }
703 
704  tmp<volScalarField> trho(this->thermo().rho());
705  const scalarField& rho = trho();
706 
707  const scalarField& T = this->thermo().T();
708  const scalarField& p = this->thermo().p();
709 
710  scalarField c0(nSpecie_);
711 
712  forAll(rho, celli)
713  {
714  scalar Ti = T[celli];
715 
716  if (Ti > Treact_)
717  {
718  const scalar rhoi = rho[celli];
719  scalar pi = p[celli];
720 
721  for (label i=0; i<nSpecie_; i++)
722  {
723  c_[i] = rhoi*Y_[i][celli]/specieThermo_[i].W();
724  c0[i] = c_[i];
725  }
726 
727  // Initialise time progress
728  scalar timeLeft = deltaT[celli];
729 
730  // Calculate the chemical source terms
731  while (timeLeft > SMALL)
732  {
733  scalar dt = timeLeft;
734  this->solve(c_, Ti, pi, dt, this->deltaTChem_[celli]);
735  timeLeft -= dt;
736  }
737 
738  deltaTMin = min(this->deltaTChem_[celli], deltaTMin);
739 
740  this->deltaTChem_[celli] =
741  min(this->deltaTChem_[celli], this->deltaTChemMax_);
742 
743  for (label i=0; i<nSpecie_; i++)
744  {
745  RR_[i][celli] =
746  (c_[i] - c0[i])*specieThermo_[i].W()/deltaT[celli];
747  }
748  }
749  else
750  {
751  for (label i=0; i<nSpecie_; i++)
752  {
753  RR_[i][celli] = 0;
754  }
755  }
756  }
757 
758  return deltaTMin;
759 }
760 
761 
762 template<class ReactionThermo, class ThermoType>
764 (
765  const scalar deltaT
766 )
767 {
768  // Don't allow the time-step to change more than a factor of 2
769  return min
770  (
772  2*deltaT
773  );
774 }
775 
776 
777 template<class ReactionThermo, class ThermoType>
779 (
780  const scalarField& deltaT
781 )
782 {
783  return this->solve<scalarField>(deltaT);
784 }
785 
786 
787 // ************************************************************************* //
Foam::constant::thermodynamic::RR
const scalar RR
Universal gas constant: default in [J/(kmol K)].
Definition: thermodynamicConstants.C:46
Foam::IOobject
Defines the attributes of an object for which implicit objectRegistry management is supported,...
Definition: IOobject.H:169
p
volScalarField & p
Definition: createFieldRefs.H:8
reactingMixture.H
s
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
Definition: gmvOutputSpray.H:25
Foam::StandardChemistryModel::~StandardChemistryModel
virtual ~StandardChemistryModel()
Destructor.
Definition: StandardChemistryModel.C:100
cp
const volScalarField & cp
Definition: setRegionSolidFields.H:8
Foam::tmp
A class for managing temporary objects.
Definition: PtrList.H:61
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::dimEnergy
const dimensionSet dimEnergy
StandardChemistryModel.H
thermo
psiReactionThermo & thermo
Definition: createFields.H:28
Foam::reactingMixture
Foam::reactingMixture.
Definition: reactingMixture.H:57
Foam::BasicChemistryModel
Basic chemistry model templated on thermodynamics.
Definition: BasicChemistryModel.H:57
thermo
Basic thermodynamics type based on the use of fitting functions for cp, h, s obtained from the templa...
Foam::endl
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:369
Foam::exp
dimensionedScalar exp(const dimensionedScalar &ds)
Definition: dimensionedScalar.C:261
Foam::StandardChemistryModel::jacobian
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &dfdc) const
Calculate the Jacobian of the system.
Definition: StandardChemistryModel.C:329
rho
rho
Definition: readInitialConditions.H:88
Foam::StandardChemistryModel::derivatives
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Calculate the derivatives in dydx.
Definition: StandardChemistryModel.C:281
composition
basicSpecieMixture & composition
Definition: createFieldRefs.H:6
Foam::min
label min(const labelHashSet &set, label minValue=labelMax)
Find the min value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:33
solve
CEqn solve()
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
trho
tmp< volScalarField > trho
Definition: setRegionSolidFields.H:4
Foam::UniformField
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:51
Foam::dimTime
const dimensionSet dimTime(0, 0, 1, 0, 0, 0, 0)
Definition: dimensionSets.H:53
R
#define R(A, B, C, D, E, F, K, M)
Foam::tmp::ref
T & ref() const
Definition: tmpI.H:227
Foam::Field< scalar >
Foam::Info
messageStream Info
Information stream (stdout output on master, null elsewhere)
correct
fvOptions correct(rho)
delta
scalar delta
Definition: LISASMDCalcMethod2.H:8
Foam::dimensionedScalar
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
Definition: dimensionedScalarFwd.H:42
Foam::StandardChemistryModel::calculate
virtual void calculate()
Calculates the reaction rates.
Definition: StandardChemistryModel.C:653
Foam::pow
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
Definition: dimensionedScalar.C:75
Foam::max
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
timeName
word timeName
Definition: getTimeIndex.H:3
Qdot
scalar Qdot
Definition: solveChemistry.H:2
Foam::StandardChemistryModel::tc
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
Definition: StandardChemistryModel.C:471
mesh
dynamicFvMesh & mesh
Definition: createDynamicFvMesh.H:6
Foam::dimMass
const dimensionSet dimMass(1, 0, 0, 0, 0, 0, 0)
Definition: dimensionSets.H:51
T
const volScalarField & T
Definition: createFieldRefs.H:2
Y
PtrList< volScalarField > & Y
Definition: createFieldRefs.H:7
Foam::StandardChemistryModel::omega
virtual void omega(const scalarField &c, const scalar T, const scalar p, scalarField &dcdt) const
dc/dt = omega, rate of change in concentration, for each species
Definition: StandardChemistryModel.C:108
Foam::SquareMatrix< scalar >
Foam::constant::mathematical::pi
constexpr scalar pi(M_PI)
Foam::ODESystem
Abstract base class for the systems of ordinary differential equations.
Definition: ODESystem.H:49
Foam::StandardChemistryModel::omegaI
virtual scalar omegaI(label iReaction, const scalarField &c, const scalar T, const scalar p, scalar &pf, scalar &cf, label &lRef, scalar &pr, scalar &cr, label &rRef) const
Return the reaction rate for iReaction and the reference.
Definition: StandardChemistryModel.C:148
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::name
word name(const expressions::valueTypeCode typeCode)
A word representation of a valueTypeCode. Empty for INVALID.
Definition: exprTraits.C:59
Foam::StandardChemistryModel::Qdot
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s3].
Definition: StandardChemistryModel.C:545
Foam::dimVolume
const dimensionSet dimVolume(pow3(dimLength))
Definition: dimensionSets.H:60
Foam::GeometricField< scalar, fvPatchField, volMesh >
Foam::Reaction
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:59
Foam::StandardChemistryModel::calculateRR
virtual tmp< volScalarField::Internal > calculateRR(const label reactionI, const label speciei) const
Return reaction rate of the speciei in reactionI.
Definition: StandardChemistryModel.C:586
UniformField.H
Foam::DimensionedField
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
Definition: DimensionedField.H:54
extrapolatedCalculatedFvPatchFields.H
Foam::StandardChemistryModel
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
Definition: StandardChemistryModel.H:61