chemistryModel.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) 2016-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 "chemistryModel.H"
27 #include "UniformField.H"
28 #include "localEulerDdtScheme.H"
29 #include "cpuLoad.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class ThermoType>
35 (
37 )
38 :
40  log_(this->lookupOrDefault("log", false)),
41  cpuLoad_(this->lookupOrDefault("cpuLoad", false)),
42  jacobianType_
43  (
44  this->found("jacobian")
45  ? jacobianTypeNames_.read(this->lookup("jacobian"))
46  : jacobianType::fast
47  ),
48  mixture_
49  (
50  dynamicCast<const multicomponentMixture<ThermoType>>(this->thermo())
51  ),
52  specieThermos_(mixture_.specieThermos()),
53  reactions_(thermo.species(), specieThermos_, this->mesh(), *this),
54  RR_(nSpecie_),
55  Y_(nSpecie_),
56  c_(nSpecie_),
57  YTpWork_(scalarField(nSpecie_ + 2)),
58  YTpYTpWork_(scalarSquareMatrix(nSpecie_ + 2)),
59  mechRedPtr_
60  (
61  chemistryReductionMethod<ThermoType>::New
62  (
63  *this,
64  *this
65  )
66  ),
67  mechRed_(*mechRedPtr_),
68  tabulationPtr_(chemistryTabulationMethod::New(*this, *this)),
69  tabulation_(*tabulationPtr_)
70 {
71  // Create the fields for the chemistry sources
72  forAll(RR_, fieldi)
73  {
74  RR_.set
75  (
76  fieldi,
78  (
79  IOobject
80  (
81  "RR." + Yvf_[fieldi].name(),
82  this->mesh().time().name(),
83  this->mesh(),
86  ),
87  thermo.mesh(),
89  )
90  );
91  }
92 
93  Info<< "chemistryModel: Number of species = " << nSpecie_
94  << " and reactions = " << nReaction() << endl;
95 
96  // When the mechanism reduction method is used, the 'active' flag for every
97  // species should be initialised (by default 'active' is true)
98  if (reduction_)
99  {
100  forAll(Yvf_, i)
101  {
103  (
104  Yvf_[i].name(),
105  this->mesh().time().name(),
106  this->mesh(),
108  );
109 
110  // Check if the species file is provided, if not set inactive
111  // and NO_WRITE
112  if (!header.headerOk())
113  {
114  this->thermo().setSpecieInactive(i);
115  }
116  }
117  }
118 
119  if (log_)
120  {
121  cpuSolveFile_ = logFile("cpu_solve.out");
122  }
123 }
124 
125 
126 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
127 
128 template<class ThermoType>
130 {}
131 
132 
133 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
134 
135 template<class ThermoType>
137 (
138  const scalar time,
139  const scalarField& YTp,
140  const label li,
141  scalarField& dYTpdt
142 ) const
143 {
144  if (reduction_)
145  {
146  forAll(sToc_, i)
147  {
148  Y_[sToc_[i]] = max(YTp[i], 0);
149  }
150  }
151  else
152  {
153  forAll(Y_, i)
154  {
155  Y_[i] = max(YTp[i], 0);
156  }
157  }
158 
159  const scalar T = YTp[nSpecie_];
160  const scalar p = YTp[nSpecie_ + 1];
161 
162  // Evaluate the mixture density
163  scalar rhoM = 0;
164  for (label i=0; i<Y_.size(); i++)
165  {
166  rhoM += Y_[i]/specieThermos_[i].rho(p, T);
167  }
168  rhoM = 1/rhoM;
169 
170  // Evaluate the concentrations
171  for (label i=0; i<Y_.size(); i ++)
172  {
173  c_[i] = rhoM/specieThermos_[i].W()*Y_[i];
174  }
175 
176  // Evaluate contributions from reactions
177  dYTpdt = Zero;
178  forAll(reactions_, ri)
179  {
180  if (!mechRed_.reactionDisabled(ri))
181  {
182  reactions_[ri].dNdtByV
183  (
184  p,
185  T,
186  c_,
187  li,
188  dYTpdt,
189  reduction_,
190  cTos_,
191  0
192  );
193  }
194  }
195 
196  // Reactions return dNdtByV, so we need to convert the result to dYdt
197  for (label i=0; i<nSpecie_; i++)
198  {
199  const scalar WiByrhoM = specieThermos_[sToc(i)].W()/rhoM;
200  scalar& dYidt = dYTpdt[i];
201  dYidt *= WiByrhoM;
202  }
203 
204  // Evaluate the effect on the thermodynamic system ...
205 
206  // Evaluate the mixture Cp
207  scalar CpM = 0;
208  for (label i=0; i<Y_.size(); i++)
209  {
210  CpM += Y_[i]*specieThermos_[i].Cp(p, T);
211  }
212 
213  // dT/dt
214  scalar& dTdt = dYTpdt[nSpecie_];
215  for (label i=0; i<nSpecie_; i++)
216  {
217  dTdt -= dYTpdt[i]*specieThermos_[sToc(i)].ha(p, T);
218  }
219  dTdt /= CpM;
220 
221  // dp/dt = 0 (pressure is assumed constant)
222  scalar& dpdt = dYTpdt[nSpecie_ + 1];
223  dpdt = 0;
224 }
225 
226 
227 template<class ThermoType>
229 (
230  const scalar t,
231  const scalarField& YTp,
232  const label li,
233  scalarField& dYTpdt,
235 ) const
236 {
237  if (reduction_)
238  {
239  forAll(sToc_, i)
240  {
241  Y_[sToc_[i]] = max(YTp[i], 0);
242  }
243  }
244  else
245  {
246  forAll(c_, i)
247  {
248  Y_[i] = max(YTp[i], 0);
249  }
250  }
251 
252  const scalar T = YTp[nSpecie_];
253  const scalar p = YTp[nSpecie_ + 1];
254 
255  // Evaluate the specific volumes and mixture density
256  scalarField& v = YTpWork_[0];
257  for (label i=0; i<Y_.size(); i++)
258  {
259  v[i] = 1/specieThermos_[i].rho(p, T);
260  }
261  scalar rhoM = 0;
262  for (label i=0; i<Y_.size(); i++)
263  {
264  rhoM += Y_[i]*v[i];
265  }
266  rhoM = 1/rhoM;
267 
268  // Evaluate the concentrations
269  for (label i=0; i<Y_.size(); i ++)
270  {
271  c_[i] = rhoM/specieThermos_[i].W()*Y_[i];
272  }
273 
274  // Evaluate the derivatives of concentration w.r.t. mass fraction
275  scalarSquareMatrix& dcdY = YTpYTpWork_[0];
276  for (label i=0; i<nSpecie_; i++)
277  {
278  const scalar rhoMByWi = rhoM/specieThermos_[sToc(i)].W();
279  switch (jacobianType_)
280  {
281  case jacobianType::fast:
282  {
283  dcdY(i, i) = rhoMByWi;
284  }
285  break;
286  case jacobianType::exact:
287  for (label j=0; j<nSpecie_; j++)
288  {
289  dcdY(i, j) =
290  rhoMByWi*((i == j) - rhoM*v[sToc(j)]*Y_[sToc(i)]);
291  }
292  break;
293  }
294  }
295 
296  // Evaluate the mixture thermal expansion coefficient
297  scalar alphavM = 0;
298  for (label i=0; i<Y_.size(); i++)
299  {
300  alphavM += Y_[i]*rhoM*v[i]*specieThermos_[i].alphav(p, T);
301  }
302 
303  // Evaluate contributions from reactions
304  dYTpdt = Zero;
305  scalarSquareMatrix& ddNdtByVdcTp = YTpYTpWork_[1];
306  for (label i=0; i<nSpecie_ + 2; i++)
307  {
308  for (label j=0; j<nSpecie_ + 2; j++)
309  {
310  ddNdtByVdcTp[i][j] = 0;
311  }
312  }
313  forAll(reactions_, ri)
314  {
315  if (!mechRed_.reactionDisabled(ri))
316  {
317  reactions_[ri].ddNdtByVdcTp
318  (
319  p,
320  T,
321  c_,
322  li,
323  dYTpdt,
324  ddNdtByVdcTp,
325  reduction_,
326  cTos_,
327  0,
328  nSpecie_,
329  YTpWork_[1],
330  YTpWork_[2]
331  );
332  }
333  }
334 
335  // Reactions return dNdtByV, so we need to convert the result to dYdt
336  for (label i=0; i<nSpecie_; i++)
337  {
338  const scalar WiByrhoM = specieThermos_[sToc(i)].W()/rhoM;
339  scalar& dYidt = dYTpdt[i];
340  dYidt *= WiByrhoM;
341 
342  for (label j=0; j<nSpecie_; j++)
343  {
344  scalar ddNidtByVdYj = 0;
345  switch (jacobianType_)
346  {
347  case jacobianType::fast:
348  {
349  const scalar ddNidtByVdcj = ddNdtByVdcTp(i, j);
350  ddNidtByVdYj = ddNidtByVdcj*dcdY(j, j);
351  }
352  break;
353  case jacobianType::exact:
354  for (label k=0; k<nSpecie_; k++)
355  {
356  const scalar ddNidtByVdck = ddNdtByVdcTp(i, k);
357  ddNidtByVdYj += ddNidtByVdck*dcdY(k, j);
358  }
359  break;
360  }
361 
362  scalar& ddYidtdYj = J(i, j);
363  ddYidtdYj = WiByrhoM*ddNidtByVdYj + rhoM*v[sToc(j)]*dYidt;
364  }
365 
366  scalar ddNidtByVdT = ddNdtByVdcTp(i, nSpecie_);
367  for (label j=0; j<nSpecie_; j++)
368  {
369  const scalar ddNidtByVdcj = ddNdtByVdcTp(i, j);
370  ddNidtByVdT -= ddNidtByVdcj*c_[sToc(j)]*alphavM;
371  }
372 
373  scalar& ddYidtdT = J(i, nSpecie_);
374  ddYidtdT = WiByrhoM*ddNidtByVdT + alphavM*dYidt;
375 
376  scalar& ddYidtdp = J(i, nSpecie_ + 1);
377  ddYidtdp = 0;
378  }
379 
380  // Evaluate the effect on the thermodynamic system ...
381 
382  // Evaluate the mixture Cp and its derivative
383  scalarField& Cp = YTpWork_[3];
384  scalar CpM = 0, dCpMdT = 0;
385  for (label i=0; i<Y_.size(); i++)
386  {
387  Cp[i] = specieThermos_[i].Cp(p, T);
388  CpM += Y_[i]*Cp[i];
389  dCpMdT += Y_[i]*specieThermos_[i].dCpdT(p, T);
390  }
391 
392  // dT/dt
393  scalarField& ha = YTpWork_[4];
394  scalar& dTdt = dYTpdt[nSpecie_];
395  for (label i=0; i<nSpecie_; i++)
396  {
397  ha[sToc(i)] = specieThermos_[sToc(i)].ha(p, T);
398  dTdt -= dYTpdt[i]*ha[sToc(i)];
399  }
400  dTdt /= CpM;
401 
402  // dp/dt = 0 (pressure is assumed constant)
403  scalar& dpdt = dYTpdt[nSpecie_ + 1];
404  dpdt = 0;
405 
406  // d(dTdt)/dY
407  for (label i=0; i<nSpecie_; i++)
408  {
409  scalar& ddTdtdYi = J(nSpecie_, i);
410  ddTdtdYi = 0;
411  for (label j=0; j<nSpecie_; j++)
412  {
413  const scalar ddYjdtdYi = J(j, i);
414  ddTdtdYi -= ddYjdtdYi*ha[sToc(j)];
415  }
416  ddTdtdYi -= Cp[sToc(i)]*dTdt;
417  ddTdtdYi /= CpM;
418  }
419 
420  // d(dTdt)/dT
421  scalar& ddTdtdT = J(nSpecie_, nSpecie_);
422  ddTdtdT = 0;
423  for (label i=0; i<nSpecie_; i++)
424  {
425  const scalar dYidt = dYTpdt[i];
426  const scalar ddYidtdT = J(i, nSpecie_);
427  ddTdtdT -= dYidt*Cp[sToc(i)] + ddYidtdT*ha[sToc(i)];
428  }
429  ddTdtdT -= dTdt*dCpMdT;
430  ddTdtdT /= CpM;
431 
432  // d(dTdt)/dp = 0 (pressure is assumed constant)
433  scalar& ddTdtdp = J(nSpecie_, nSpecie_ + 1);
434  ddTdtdp = 0;
435 
436  // d(dpdt)/dYiTp = 0 (pressure is assumed constant)
437  for (label i=0; i<nSpecie_ + 2; i++)
438  {
439  scalar& ddpdtdYiTp = J(nSpecie_ + 1, i);
440  ddpdtdYiTp = 0;
441  }
442 }
443 
444 
445 template<class ThermoType>
448 (
449  const label reactioni
450 ) const
451 {
453  for (label i=0; i<nSpecie_; i++)
454  {
455  RR.set
456  (
457  i,
459  (
460  "RR." + Yvf_[i].name(),
461  this->mesh(),
463  ).ptr()
464  );
465  }
466 
467  if (!this->chemistry_ || mechRed_.reactionDisabled(reactioni))
468  {
469  return RR;
470  }
471 
472  tmp<volScalarField> trhovf(this->thermo().rho());
473  const volScalarField& rhovf = trhovf();
474 
475  const volScalarField& Tvf = this->thermo().T();
476  const volScalarField& pvf = this->thermo().p();
477 
478  scalarField& dNdtByV = YTpWork_[0];
479 
480  reactionEvaluationScope scope(*this);
481 
482  const Reaction<ThermoType>& R = reactions_[reactioni];
483 
484  forAll(rhovf, celli)
485  {
486  const scalar rho = rhovf[celli];
487  const scalar T = Tvf[celli];
488  const scalar p = pvf[celli];
489 
490  for (label i=0; i<nSpecie_; i++)
491  {
492  const scalar Yi = Yvf_[i][celli];
493  c_[i] = rho*Yi/specieThermos_[i].W();
494  }
495 
496  dNdtByV = Zero;
497 
498  R.dNdtByV
499  (
500  p,
501  T,
502  c_,
503  celli,
504  dNdtByV,
505  reduction_,
506  cTos_,
507  0
508  );
509 
510  for (label i=0; i<nSpecie_; i++)
511  {
512  RR[i][celli] = dNdtByV[i]*specieThermos_[i].W();
513  }
514  }
515 
516  return RR;
517 }
518 
519 
520 template<class ThermoType>
522 {
523  if (!this->chemistry_)
524  {
525  return;
526  }
527 
528  tmp<volScalarField> trhovf(this->thermo().rho());
529  const volScalarField& rhovf = trhovf();
530 
531  const volScalarField& Tvf = this->thermo().T();
532  const volScalarField& pvf = this->thermo().p();
533 
534  scalarField& dNdtByV = YTpWork_[0];
535 
536  reactionEvaluationScope scope(*this);
537 
538  forAll(rhovf, celli)
539  {
540  const scalar rho = rhovf[celli];
541  const scalar T = Tvf[celli];
542  const scalar p = pvf[celli];
543 
544  for (label i=0; i<nSpecie_; i++)
545  {
546  const scalar Yi = Yvf_[i][celli];
547  c_[i] = rho*Yi/specieThermos_[i].W();
548  }
549 
550  dNdtByV = Zero;
551 
552  forAll(reactions_, ri)
553  {
554  if (!mechRed_.reactionDisabled(ri))
555  {
556  reactions_[ri].dNdtByV
557  (
558  p,
559  T,
560  c_,
561  celli,
562  dNdtByV,
563  reduction_,
564  cTos_,
565  0
566  );
567  }
568  }
569 
570  for (label i=0; i<mechRed_.nActiveSpecies(); i++)
571  {
572  RR_[sToc(i)][celli] = dNdtByV[i]*specieThermos_[sToc(i)].W();
573  }
574  }
575 }
576 
577 
578 template<class ThermoType>
579 template<class DeltaTType>
581 (
582  const DeltaTType& deltaT
583 )
584 {
585  optionalCpuLoad& chemistryCpuLoad
586  (
587  optionalCpuLoad::New(name() + ":cpuLoad", this->mesh(), cpuLoad_)
588  );
589 
590  // CPU time logging
591  cpuTime solveCpuTime;
592  scalar totalSolveCpuTime = 0;
593 
594  if (!this->chemistry_)
595  {
596  return great;
597  }
598 
599  const volScalarField& rho0vf =
600  this->mesh().template lookupObject<volScalarField>
601  (
602  this->thermo().phasePropertyName("rho")
603  ).oldTime();
604 
605  const volScalarField& T0vf = this->thermo().T().oldTime();
606  const volScalarField& p0vf = this->thermo().p().oldTime();
607 
608  reactionEvaluationScope scope(*this);
609 
610  scalarField Y0(nSpecie_);
611 
612  // Composition vector (Yi, T, p, deltaT)
613  scalarField phiq(nEqns() + 1);
614  scalarField Rphiq(nEqns() + 1);
615 
616  // Minimum chemical timestep
617  scalar deltaTMin = great;
618 
619  tabulation_.reset();
620  chemistryCpuLoad.resetCpuTime();
621 
622  forAll(rho0vf, celli)
623  {
624  const scalar rho0 = rho0vf[celli];
625 
626  scalar p = p0vf[celli];
627  scalar T = T0vf[celli];
628 
629  for (label i=0; i<nSpecie_; i++)
630  {
631  Y_[i] = Y0[i] = Yvf_[i].oldTime()[celli];
632  }
633 
634  for (label i=0; i<nSpecie_; i++)
635  {
636  phiq[i] = Yvf_[i].oldTime()[celli];
637  }
638  phiq[nSpecie()] = T;
639  phiq[nSpecie() + 1] = p;
640  phiq[nSpecie() + 2] = deltaT[celli];
641 
642  // Initialise time progress
643  scalar timeLeft = deltaT[celli];
644 
645  // Not sure if this is necessary
646  Rphiq = Zero;
647 
648  // When tabulation is active (short-circuit evaluation for retrieve)
649  // It first tries to retrieve the solution of the system with the
650  // information stored through the tabulation method
651  if (tabulation_.retrieve(phiq, Rphiq))
652  {
653  // Retrieved solution stored in Rphiq
654  for (label i=0; i<nSpecie(); i++)
655  {
656  Y_[i] = Rphiq[i];
657  }
658  T = Rphiq[nSpecie()];
659  p = Rphiq[nSpecie() + 1];
660  }
661  // This position is reached when tabulation is not used OR
662  // if the solution is not retrieved.
663  // In the latter case, it adds the information to the tabulation
664  // (it will either expand the current data or add a new stored point).
665  else
666  {
667  if (reduction_)
668  {
669  // Compute concentrations
670  for (label i=0; i<nSpecie_; i++)
671  {
672  c_[i] = rho0*Y_[i]/specieThermos_[i].W();
673  }
674 
675  // Reduce mechanism change the number of species (only active)
676  mechRed_.reduceMechanism(p, T, c_, cTos_, sToc_, celli);
677 
678  // Set the simplified mass fraction field
679  sY_.setSize(nSpecie_);
680  for (label i=0; i<nSpecie_; i++)
681  {
682  sY_[i] = Y_[sToc(i)];
683  }
684  }
685 
686  if (log_)
687  {
688  // Reset the solve time
689  solveCpuTime.cpuTimeIncrement();
690  }
691 
692  // Calculate the chemical source terms
693  while (timeLeft > small)
694  {
695  scalar dt = timeLeft;
696  if (reduction_)
697  {
698  // Solve the reduced set of ODE
699  solve
700  (
701  p,
702  T,
703  sY_,
704  celli,
705  dt,
706  deltaTChem_[celli]
707  );
708 
709  for (label i=0; i<mechRed_.nActiveSpecies(); i++)
710  {
711  Y_[sToc_[i]] = sY_[i];
712  }
713  }
714  else
715  {
716  solve(p, T, Y_, celli, dt, deltaTChem_[celli]);
717  }
718  timeLeft -= dt;
719  }
720 
721  if (log_)
722  {
723  totalSolveCpuTime += solveCpuTime.cpuTimeIncrement();
724  }
725 
726  // If tabulation is used, we add the information computed here to
727  // the stored points (either expand or add)
728  if (tabulation_.tabulates())
729  {
730  forAll(Y_, i)
731  {
732  Rphiq[i] = Y_[i];
733  }
734  Rphiq[Rphiq.size()-3] = T;
735  Rphiq[Rphiq.size()-2] = p;
736  Rphiq[Rphiq.size()-1] = deltaT[celli];
737 
738  tabulation_.add
739  (
740  phiq,
741  Rphiq,
742  mechRed_.nActiveSpecies(),
743  celli,
744  deltaT[celli]
745  );
746  }
747 
748  // When operations are done and if mechanism reduction is active,
749  // the number of species (which also affects nEqns) is set back
750  // to the total number of species (stored in the mechRed object)
751  if (reduction_)
752  {
753  setNSpecie(mechRed_.nSpecie());
754  }
755 
756  deltaTMin = min(deltaTChem_[celli], deltaTMin);
757  deltaTChem_[celli] = min(deltaTChem_[celli], deltaTChemMax_);
758  }
759 
760  // Set the RR vector (used in the solver)
761  for (label i=0; i<nSpecie_; i++)
762  {
763  RR_[i][celli] = rho0*(Y_[i] - Y0[i])/deltaT[celli];
764  }
765 
766  if (cpuLoad_)
767  {
768  chemistryCpuLoad.cpuTimeIncrement(celli);
769  }
770  }
771 
772  if (log_)
773  {
774  cpuSolveFile_()
775  << this->time().userTimeValue()
776  << " " << totalSolveCpuTime << endl;
777  }
778 
779  mechRed_.update();
780  tabulation_.update();
781 
782  if (reduction_)
783  {
784  this->thermo().syncSpeciesActive();
785  }
786 
787  return deltaTMin;
788 }
789 
790 
791 template<class ThermoType>
793 (
794  const scalar deltaT
795 )
796 {
797  // Don't allow the time-step to change more than a factor of 2
798  return min
799  (
801  2*deltaT
802  );
803 }
804 
805 
806 template<class ThermoType>
808 (
809  const scalarField& deltaT
810 )
811 {
812  return this->solve<scalarField>(deltaT);
813 }
814 
815 
816 template<class ThermoType>
819 {
821  (
823  (
824  "tc",
825  this->mesh(),
826  dimensionedScalar(dimTime, small),
827  extrapolatedCalculatedFvPatchScalarField::typeName
828  )
829  );
830  scalarField& tc = ttc.ref();
831 
832  if (!this->chemistry_)
833  {
834  ttc.ref().correctBoundaryConditions();
835  return ttc;
836  }
837 
838  tmp<volScalarField> trhovf(this->thermo().rho());
839  const volScalarField& rhovf = trhovf();
840 
841  const volScalarField& Tvf = this->thermo().T();
842  const volScalarField& pvf = this->thermo().p();
843 
844  reactionEvaluationScope scope(*this);
845 
846  forAll(rhovf, celli)
847  {
848  const scalar rho = rhovf[celli];
849  const scalar T = Tvf[celli];
850  const scalar p = pvf[celli];
851 
852  for (label i=0; i<nSpecie_; i++)
853  {
854  c_[i] = rho*Yvf_[i][celli]/specieThermos_[i].W();
855  }
856 
857  // A reaction's rate scale is calculated as its molar
858  // production rate divided by the total number of moles in the
859  // system.
860  //
861  // The system rate scale is the average of the reactions' rate
862  // scales weighted by the reactions' molar production rates. This
863  // weighting ensures that dominant reactions provide the largest
864  // contribution to the system rate scale.
865  //
866  // The system time scale is then the reciprocal of the system rate
867  // scale.
868  //
869  // Contributions from forward and reverse reaction rates are
870  // handled independently and identically so that reversible
871  // reactions produce the same result as the equivalent pair of
872  // irreversible reactions.
873 
874  scalar sumW = 0, sumWRateByCTot = 0;
875  forAll(reactions_, i)
876  {
877  const Reaction<ThermoType>& R = reactions_[i];
878  scalar omegaf, omegar;
879  R.omega(p, T, c_, celli, omegaf, omegar);
880 
881  scalar wf = 0;
882  forAll(R.rhs(), s)
883  {
884  wf += R.rhs()[s].stoichCoeff*omegaf;
885  }
886  sumW += wf;
887  sumWRateByCTot += sqr(wf);
888 
889  scalar wr = 0;
890  forAll(R.lhs(), s)
891  {
892  wr += R.lhs()[s].stoichCoeff*omegar;
893  }
894  sumW += wr;
895  sumWRateByCTot += sqr(wr);
896  }
897 
898  tc[celli] =
899  sumWRateByCTot == 0 ? vGreat : sumW/sumWRateByCTot*sum(c_);
900  }
901 
902  ttc.ref().correctBoundaryConditions();
903  return ttc;
904 }
905 
906 
907 template<class ThermoType>
910 {
911  tmp<volScalarField> tQdot
912  (
914  (
915  "Qdot",
916  this->mesh_,
918  )
919  );
920 
921  if (!this->chemistry_)
922  {
923  return tQdot;
924  }
925 
926  reactionEvaluationScope scope(*this);
927 
928  scalarField& Qdot = tQdot.ref();
929 
930  forAll(Yvf_, i)
931  {
932  forAll(Qdot, celli)
933  {
934  const scalar hi = specieThermos_[i].hf();
935  Qdot[celli] -= hi*RR_[i][celli];
936  }
937  }
938 
939  return tQdot;
940 }
941 
942 
943 // ************************************************************************* //
scalar Cp(const scalar p, const scalar T) const
Definition: EtoHthermo.H:2
scalar ha(const scalar p, const scalar T) const
Definition: EtoHthermo.H:20
label k
bool found
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh > > New(const word &name, const Mesh &mesh, const dimensionSet &, const Field< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
static tmp< GeometricField< Type, PatchField, GeoMesh > > New(const word &name, const Internal &, const PtrList< PatchField< Type >> &, const HashPtrTable< Source > &=HashPtrTable< Source >())
Return a temporary field constructed from name,.
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:99
const Time & time() const
Return time.
Definition: IOobject.C:318
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const FieldType & oldTime() const
Return the old time field.
Definition: OldTimeField.C:315
A templated 1D list of pointers to objects of type <T>, where the size of the array is known and used...
Definition: PtrList.H:75
Simple extension of ThermoType to handle reaction kinetics in addition to the equilibrium thermodynam...
Definition: Reaction.H:72
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:48
const fluidMulticomponentThermo & thermo() const
Return const access to the thermo.
jacobianType
Enumeration for the type of Jacobian to be calculated by the.
const fvMesh & mesh() const
Return const access to the mesh.
virtual const fvMesh & mesh() const =0
Return const access to the mesh.
virtual const volScalarField & T() const =0
Temperature [K].
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
virtual void jacobian(const scalar t, const scalarField &YTp, const label li, scalarField &dYTpdt, scalarSquareMatrix &J) const
Calculate the ODE jacobian.
virtual tmp< volScalarField > tc() const
Return the chemical time scale.
virtual void derivatives(const scalar t, const scalarField &YTp, const label li, scalarField &dYTpdt) const
Calculate the ODE derivatives.
virtual label nReaction() const
The number of reactions.
virtual ~chemistryModel()
Destructor.
virtual PtrList< volScalarField::Internal > reactionRR(const label reactioni) const
Return reaction rates of the species in reactioni [kg/m^3/s].
virtual tmp< volScalarField > Qdot() const
Return the heat release rate [kg/m/s^3].
virtual void calculate()
Calculates the reaction rates.
An abstract class for methods of chemical mechanism reduction.
An abstract class for chemistry tabulation.
Starts timing CPU usage and return elapsed time from start.
Definition: cpuTime.H:55
double cpuTimeIncrement() const
Return CPU time (in seconds) since last call to cpuTimeIncrement()
Definition: cpuTime.C:60
const fileName & name() const
Return the dictionary name.
Definition: dictionary.H:111
Base-class for multi-component fluid thermodynamic properties.
virtual const volScalarField & p() const =0
Pressure [Pa].
Foam::multicomponentMixture.
void syncSpeciesActive() const
Synchronise the specie active flags.
void setSpecieInactive(const label speciei) const
Set specie inactive.
Extends base chemistry model adding an ODESystem and the reduction maps needed for tabulation.
autoPtr< OFstream > logFile(const word &name) const
Create and return a TDAC log file of the given name.
virtual void resetCpuTime()
Reset the CPU time (dummy)
Definition: cpuLoad.H:101
virtual void cpuTimeIncrement(const label celli)
Cache the CPU time increment for celli (dummy)
Definition: cpuLoad.H:105
static optionalCpuLoad & New(const word &name, const polyMesh &mesh, const bool loadBalancing)
Construct from polyMesh if loadBalancing is true.
Definition: cpuLoad.C:63
A class for managing temporary objects.
Definition: tmp.H:55
T & ref() const
Return non-const reference or generate a fatal error.
Definition: tmpI.H:181
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:531
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
void read(Istream &, label &, const dictionary &)
In-place read with dictionary lookup.
autoPtr< CompressibleMomentumTransportModel > New(const volScalarField &rho, const volVectorField &U, const surfaceScalarField &phi, const viscosity &viscosity)
const dimensionedScalar RR
Universal gas constant: default SI units: [J/kmol/K].
static const zero Zero
Definition: zero.H:97
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
const dimensionSet dimEnergy
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
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh > &df)
messageStream Info
To & dynamicCast(From &r)
Reference type cast template function,.
Definition: typeInfo.H:91
layerAndWeight min(const layerAndWeight &a, const layerAndWeight &b)
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
const dimensionSet dimMass
void T(FieldField< Field, Type > &f1, const FieldField< Field, Type > &f2)
dimensioned< scalar > dimensionedScalar
Dimensioned scalar obtained from generic dimensioned type.
SolverPerformance< Type > solve(fvMatrix< Type > &, const word &)
Solve returning the solution statistics given convergence tolerance.
scalar rho0
const label nSpecie
scalarList Y0(nSpecie, 0.0)
scalar Qdot
Definition: solveChemistry.H:2
volScalarField & p
fluidMulticomponentThermo & thermo
Definition: createFields.H:31