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-2025 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  this->thermo().syncSpeciesActive();
119  }
120 
121  if (log_)
122  {
123  cpuSolveFile_ = logFile("cpu_solve.out");
124  }
125 }
126 
127 
128 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
129 
130 template<class ThermoType>
132 {}
133 
134 
135 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
136 
137 template<class ThermoType>
139 (
140  const scalar time,
141  const scalarField& YTp,
142  const label li,
143  scalarField& dYTpdt
144 ) const
145 {
146  if (reduction_)
147  {
148  forAll(sToc_, i)
149  {
150  Y_[sToc_[i]] = max(YTp[i], 0);
151  }
152  }
153  else
154  {
155  forAll(Y_, i)
156  {
157  Y_[i] = max(YTp[i], 0);
158  }
159  }
160 
161  const scalar T = YTp[nSpecie_];
162  const scalar p = YTp[nSpecie_ + 1];
163 
164  // Evaluate the mixture density
165  scalar rhoM = 0;
166  for (label i=0; i<Y_.size(); i++)
167  {
168  rhoM += Y_[i]/specieThermos_[i].rho(p, T);
169  }
170  rhoM = 1/rhoM;
171 
172  // Evaluate the concentrations
173  for (label i=0; i<Y_.size(); i ++)
174  {
175  c_[i] = rhoM/specieThermos_[i].W()*Y_[i];
176  }
177 
178  // Evaluate contributions from reactions
179  dYTpdt = Zero;
180  forAll(reactions_, ri)
181  {
182  if (!mechRed_.reactionDisabled(ri))
183  {
184  reactions_[ri].dNdtByV
185  (
186  p,
187  T,
188  c_,
189  li,
190  dYTpdt,
191  reduction_,
192  cTos_,
193  0
194  );
195  }
196  }
197 
198  // Reactions return dNdtByV, so we need to convert the result to dYdt
199  for (label i=0; i<nSpecie_; i++)
200  {
201  const scalar WiByrhoM = specieThermos_[sToc(i)].W()/rhoM;
202  scalar& dYidt = dYTpdt[i];
203  dYidt *= WiByrhoM;
204  }
205 
206  // Evaluate the effect on the thermodynamic system ...
207 
208  // Evaluate the mixture Cp
209  scalar CpM = 0;
210  for (label i=0; i<Y_.size(); i++)
211  {
212  CpM += Y_[i]*specieThermos_[i].Cp(p, T);
213  }
214 
215  // dT/dt
216  scalar& dTdt = dYTpdt[nSpecie_];
217  for (label i=0; i<nSpecie_; i++)
218  {
219  dTdt -= dYTpdt[i]*specieThermos_[sToc(i)].ha(p, T);
220  }
221  dTdt /= CpM;
222 
223  // dp/dt = 0 (pressure is assumed constant)
224  scalar& dpdt = dYTpdt[nSpecie_ + 1];
225  dpdt = 0;
226 }
227 
228 
229 template<class ThermoType>
231 (
232  const scalar t,
233  const scalarField& YTp,
234  const label li,
235  scalarField& dYTpdt,
237 ) const
238 {
239  if (reduction_)
240  {
241  forAll(sToc_, i)
242  {
243  Y_[sToc_[i]] = max(YTp[i], 0);
244  }
245  }
246  else
247  {
248  forAll(c_, i)
249  {
250  Y_[i] = max(YTp[i], 0);
251  }
252  }
253 
254  const scalar T = YTp[nSpecie_];
255  const scalar p = YTp[nSpecie_ + 1];
256 
257  // Evaluate the specific volumes and mixture density
258  scalarField& v = YTpWork_[0];
259  for (label i=0; i<Y_.size(); i++)
260  {
261  v[i] = 1/specieThermos_[i].rho(p, T);
262  }
263  scalar rhoM = 0;
264  for (label i=0; i<Y_.size(); i++)
265  {
266  rhoM += Y_[i]*v[i];
267  }
268  rhoM = 1/rhoM;
269 
270  // Evaluate the concentrations
271  for (label i=0; i<Y_.size(); i ++)
272  {
273  c_[i] = rhoM/specieThermos_[i].W()*Y_[i];
274  }
275 
276  // Evaluate the derivatives of concentration w.r.t. mass fraction
277  scalarSquareMatrix& dcdY = YTpYTpWork_[0];
278  for (label i=0; i<nSpecie_; i++)
279  {
280  const scalar rhoMByWi = rhoM/specieThermos_[sToc(i)].W();
281  switch (jacobianType_)
282  {
283  case jacobianType::fast:
284  {
285  dcdY(i, i) = rhoMByWi;
286  }
287  break;
288  case jacobianType::exact:
289  for (label j=0; j<nSpecie_; j++)
290  {
291  dcdY(i, j) =
292  rhoMByWi*((i == j) - rhoM*v[sToc(j)]*Y_[sToc(i)]);
293  }
294  break;
295  }
296  }
297 
298  // Evaluate the mixture thermal expansion coefficient
299  scalar alphavM = 0;
300  for (label i=0; i<Y_.size(); i++)
301  {
302  alphavM += Y_[i]*rhoM*v[i]*specieThermos_[i].alphav(p, T);
303  }
304 
305  // Evaluate contributions from reactions
306  dYTpdt = Zero;
307  scalarSquareMatrix& ddNdtByVdcTp = YTpYTpWork_[1];
308  for (label i=0; i<nSpecie_ + 2; i++)
309  {
310  for (label j=0; j<nSpecie_ + 2; j++)
311  {
312  ddNdtByVdcTp[i][j] = 0;
313  }
314  }
315  forAll(reactions_, ri)
316  {
317  if (!mechRed_.reactionDisabled(ri))
318  {
319  reactions_[ri].ddNdtByVdcTp
320  (
321  p,
322  T,
323  c_,
324  li,
325  dYTpdt,
326  ddNdtByVdcTp,
327  reduction_,
328  cTos_,
329  0,
330  nSpecie_,
331  YTpWork_[1],
332  YTpWork_[2]
333  );
334  }
335  }
336 
337  // Reactions return dNdtByV, so we need to convert the result to dYdt
338  for (label i=0; i<nSpecie_; i++)
339  {
340  const scalar WiByrhoM = specieThermos_[sToc(i)].W()/rhoM;
341  scalar& dYidt = dYTpdt[i];
342  dYidt *= WiByrhoM;
343 
344  for (label j=0; j<nSpecie_; j++)
345  {
346  scalar ddNidtByVdYj = 0;
347  switch (jacobianType_)
348  {
349  case jacobianType::fast:
350  {
351  const scalar ddNidtByVdcj = ddNdtByVdcTp(i, j);
352  ddNidtByVdYj = ddNidtByVdcj*dcdY(j, j);
353  }
354  break;
355  case jacobianType::exact:
356  for (label k=0; k<nSpecie_; k++)
357  {
358  const scalar ddNidtByVdck = ddNdtByVdcTp(i, k);
359  ddNidtByVdYj += ddNidtByVdck*dcdY(k, j);
360  }
361  break;
362  }
363 
364  scalar& ddYidtdYj = J(i, j);
365  ddYidtdYj = WiByrhoM*ddNidtByVdYj + rhoM*v[sToc(j)]*dYidt;
366  }
367 
368  scalar ddNidtByVdT = ddNdtByVdcTp(i, nSpecie_);
369  for (label j=0; j<nSpecie_; j++)
370  {
371  const scalar ddNidtByVdcj = ddNdtByVdcTp(i, j);
372  ddNidtByVdT -= ddNidtByVdcj*c_[sToc(j)]*alphavM;
373  }
374 
375  scalar& ddYidtdT = J(i, nSpecie_);
376  ddYidtdT = WiByrhoM*ddNidtByVdT + alphavM*dYidt;
377 
378  scalar& ddYidtdp = J(i, nSpecie_ + 1);
379  ddYidtdp = 0;
380  }
381 
382  // Evaluate the effect on the thermodynamic system ...
383 
384  // Evaluate the mixture Cp and its derivative
385  scalarField& Cp = YTpWork_[3];
386  scalar CpM = 0, dCpMdT = 0;
387  for (label i=0; i<Y_.size(); i++)
388  {
389  Cp[i] = specieThermos_[i].Cp(p, T);
390  CpM += Y_[i]*Cp[i];
391  dCpMdT += Y_[i]*specieThermos_[i].dCpdT(p, T);
392  }
393 
394  // dT/dt
395  scalarField& ha = YTpWork_[4];
396  scalar& dTdt = dYTpdt[nSpecie_];
397  for (label i=0; i<nSpecie_; i++)
398  {
399  ha[sToc(i)] = specieThermos_[sToc(i)].ha(p, T);
400  dTdt -= dYTpdt[i]*ha[sToc(i)];
401  }
402  dTdt /= CpM;
403 
404  // dp/dt = 0 (pressure is assumed constant)
405  scalar& dpdt = dYTpdt[nSpecie_ + 1];
406  dpdt = 0;
407 
408  // d(dTdt)/dY
409  for (label i=0; i<nSpecie_; i++)
410  {
411  scalar& ddTdtdYi = J(nSpecie_, i);
412  ddTdtdYi = 0;
413  for (label j=0; j<nSpecie_; j++)
414  {
415  const scalar ddYjdtdYi = J(j, i);
416  ddTdtdYi -= ddYjdtdYi*ha[sToc(j)];
417  }
418  ddTdtdYi -= Cp[sToc(i)]*dTdt;
419  ddTdtdYi /= CpM;
420  }
421 
422  // d(dTdt)/dT
423  scalar& ddTdtdT = J(nSpecie_, nSpecie_);
424  ddTdtdT = 0;
425  for (label i=0; i<nSpecie_; i++)
426  {
427  const scalar dYidt = dYTpdt[i];
428  const scalar ddYidtdT = J(i, nSpecie_);
429  ddTdtdT -= dYidt*Cp[sToc(i)] + ddYidtdT*ha[sToc(i)];
430  }
431  ddTdtdT -= dTdt*dCpMdT;
432  ddTdtdT /= CpM;
433 
434  // d(dTdt)/dp = 0 (pressure is assumed constant)
435  scalar& ddTdtdp = J(nSpecie_, nSpecie_ + 1);
436  ddTdtdp = 0;
437 
438  // d(dpdt)/dYiTp = 0 (pressure is assumed constant)
439  for (label i=0; i<nSpecie_ + 2; i++)
440  {
441  scalar& ddpdtdYiTp = J(nSpecie_ + 1, i);
442  ddpdtdYiTp = 0;
443  }
444 }
445 
446 
447 template<class ThermoType>
450 (
451  const label reactioni
452 ) const
453 {
456  (
457  "RR:" + reactions_[reactioni].name(),
458  this->mesh(),
460  );
462 
463  if (!this->chemistry_ || mechRed_.reactionDisabled(reactioni))
464  {
465  return tRR;
466  }
467 
468  tmp<volScalarField> trhovf(this->thermo().rho());
469  const volScalarField& rhovf = trhovf();
470 
471  const volScalarField& Tvf = this->thermo().T();
472  const volScalarField& pvf = this->thermo().p();
473 
474  reactionEvaluationScope scope(*this);
475 
476  const Reaction<ThermoType>& R = reactions_[reactioni];
477 
478  forAll(rhovf, celli)
479  {
480  const scalar rho = rhovf[celli];
481  const scalar T = Tvf[celli];
482  const scalar p = pvf[celli];
483 
484  for (label i=0; i<nSpecie_; i++)
485  {
486  const scalar Yi = Yvf_[i][celli];
487  c_[i] = rho*Yi/specieThermos_[i].W();
488  }
489 
490  scalar omegaf, omegar;
491 
492  RR[celli] =
493  R.omega
494  (
495  p,
496  T,
497  c_,
498  celli,
499  omegaf,
500  omegar
501  );
502  }
503 
504  return tRR;
505 }
506 
507 
508 template<class ThermoType>
511 (
512  const label reactioni
513 ) const
514 {
516  for (label i=0; i<nSpecie_; i++)
517  {
518  RR.set
519  (
520  i,
522  (
523  "RR:" + reactions_[reactioni].name() + ":" + Yvf_[i].name(),
524  this->mesh(),
526  ).ptr()
527  );
528  }
529 
530  if (!this->chemistry_ || mechRed_.reactionDisabled(reactioni))
531  {
532  return RR;
533  }
534 
535  tmp<volScalarField> trhovf(this->thermo().rho());
536  const volScalarField& rhovf = trhovf();
537 
538  const volScalarField& Tvf = this->thermo().T();
539  const volScalarField& pvf = this->thermo().p();
540 
541  scalarField& dNdtByV = YTpWork_[0];
542 
543  reactionEvaluationScope scope(*this);
544 
545  const Reaction<ThermoType>& R = reactions_[reactioni];
546 
547  forAll(rhovf, celli)
548  {
549  const scalar rho = rhovf[celli];
550  const scalar T = Tvf[celli];
551  const scalar p = pvf[celli];
552 
553  for (label i=0; i<nSpecie_; i++)
554  {
555  const scalar Yi = Yvf_[i][celli];
556  c_[i] = rho*Yi/specieThermos_[i].W();
557  }
558 
559  dNdtByV = Zero;
560 
561  R.dNdtByV
562  (
563  p,
564  T,
565  c_,
566  celli,
567  dNdtByV,
568  reduction_,
569  cTos_,
570  0
571  );
572 
573  for (label i=0; i<nSpecie_; i++)
574  {
575  RR[i][celli] = dNdtByV[i]*specieThermos_[i].W();
576  }
577  }
578 
579  return RR;
580 }
581 
582 
583 template<class ThermoType>
585 {
586  if (!this->chemistry_)
587  {
588  return;
589  }
590 
591  tmp<volScalarField> trhovf(this->thermo().rho());
592  const volScalarField& rhovf = trhovf();
593 
594  const volScalarField& Tvf = this->thermo().T();
595  const volScalarField& pvf = this->thermo().p();
596 
597  scalarField& dNdtByV = YTpWork_[0];
598 
599  reactionEvaluationScope scope(*this);
600 
601  forAll(rhovf, celli)
602  {
603  const scalar rho = rhovf[celli];
604  const scalar T = Tvf[celli];
605  const scalar p = pvf[celli];
606 
607  for (label i=0; i<nSpecie_; i++)
608  {
609  const scalar Yi = Yvf_[i][celli];
610  c_[i] = rho*Yi/specieThermos_[i].W();
611  }
612 
613  dNdtByV = Zero;
614 
615  forAll(reactions_, ri)
616  {
617  if (!mechRed_.reactionDisabled(ri))
618  {
619  reactions_[ri].dNdtByV
620  (
621  p,
622  T,
623  c_,
624  celli,
625  dNdtByV,
626  reduction_,
627  cTos_,
628  0
629  );
630  }
631  }
632 
633  for (label i=0; i<mechRed_.nActiveSpecies(); i++)
634  {
635  RR_[sToc(i)][celli] = dNdtByV[i]*specieThermos_[sToc(i)].W();
636  }
637  }
638 }
639 
640 
641 template<class ThermoType>
642 template<class DeltaTType>
644 (
645  const DeltaTType& deltaT
646 )
647 {
648  optionalCpuLoad& chemistryCpuLoad
649  (
650  optionalCpuLoad::New(name() + ":cpuLoad", this->mesh(), cpuLoad_)
651  );
652 
653  // CPU time logging
654  cpuTime solveCpuTime;
655  scalar totalSolveCpuTime = 0;
656 
657  if (!this->chemistry_)
658  {
659  return great;
660  }
661 
662  const volScalarField& rho0vf =
663  this->mesh().template lookupObject<volScalarField>
664  (
665  this->thermo().phasePropertyName("rho")
666  ).oldTime();
667 
668  const volScalarField& T0vf = this->thermo().T().oldTime();
669  const volScalarField& p0vf = this->thermo().p().oldTime();
670 
671  reactionEvaluationScope scope(*this);
672 
673  scalarField Y0(nSpecie_);
674 
675  // Composition vector (Yi, T, p, deltaT)
676  scalarField phiq(nEqns() + 1);
677  scalarField Rphiq(nEqns() + 1);
678 
679  // Minimum chemical timestep
680  scalar deltaTMin = great;
681 
682  tabulation_.reset();
683  chemistryCpuLoad.resetCpuTime();
684 
685  forAll(rho0vf, celli)
686  {
687  const scalar rho0 = rho0vf[celli];
688 
689  scalar p = p0vf[celli];
690  scalar T = T0vf[celli];
691 
692  for (label i=0; i<nSpecie_; i++)
693  {
694  Y_[i] = Y0[i] = Yvf_[i].oldTime()[celli];
695  }
696 
697  for (label i=0; i<nSpecie_; i++)
698  {
699  phiq[i] = Yvf_[i].oldTime()[celli];
700  }
701  phiq[nSpecie()] = T;
702  phiq[nSpecie() + 1] = p;
703  phiq[nSpecie() + 2] = deltaT[celli];
704 
705  // Initialise time progress
706  scalar timeLeft = deltaT[celli];
707 
708  // Not sure if this is necessary
709  Rphiq = Zero;
710 
711  // When tabulation is active (short-circuit evaluation for retrieve)
712  // It first tries to retrieve the solution of the system with the
713  // information stored through the tabulation method
714  if (tabulation_.retrieve(phiq, Rphiq))
715  {
716  // Retrieved solution stored in Rphiq
717  for (label i=0; i<nSpecie(); i++)
718  {
719  Y_[i] = Rphiq[i];
720  }
721  T = Rphiq[nSpecie()];
722  p = Rphiq[nSpecie() + 1];
723  }
724  // This position is reached when tabulation is not used OR
725  // if the solution is not retrieved.
726  // In the latter case, it adds the information to the tabulation
727  // (it will either expand the current data or add a new stored point).
728  else
729  {
730  if (reduction_)
731  {
732  // Compute concentrations
733  for (label i=0; i<nSpecie_; i++)
734  {
735  c_[i] = rho0*Y_[i]/specieThermos_[i].W();
736  }
737 
738  // Reduce mechanism change the number of species (only active)
739  mechRed_.reduceMechanism(p, T, c_, cTos_, sToc_, celli);
740 
741  // Set the simplified mass fraction field
742  sY_.setSize(nSpecie_);
743  for (label i=0; i<nSpecie_; i++)
744  {
745  sY_[i] = Y_[sToc(i)];
746  }
747  }
748 
749  if (log_)
750  {
751  // Reset the solve time
752  solveCpuTime.cpuTimeIncrement();
753  }
754 
755  // Calculate the chemical source terms
756  while (timeLeft > small)
757  {
758  scalar dt = timeLeft;
759  if (reduction_)
760  {
761  // Solve the reduced set of ODE
762  solve
763  (
764  p,
765  T,
766  sY_,
767  celli,
768  dt,
769  deltaTChem_[celli]
770  );
771 
772  for (label i=0; i<mechRed_.nActiveSpecies(); i++)
773  {
774  Y_[sToc_[i]] = sY_[i];
775  }
776  }
777  else
778  {
779  solve(p, T, Y_, celli, dt, deltaTChem_[celli]);
780  }
781  timeLeft -= dt;
782  }
783 
784  if (log_)
785  {
786  totalSolveCpuTime += solveCpuTime.cpuTimeIncrement();
787  }
788 
789  // If tabulation is used, we add the information computed here to
790  // the stored points (either expand or add)
791  if (tabulation_.tabulates())
792  {
793  forAll(Y_, i)
794  {
795  Rphiq[i] = Y_[i];
796  }
797  Rphiq[Rphiq.size()-3] = T;
798  Rphiq[Rphiq.size()-2] = p;
799  Rphiq[Rphiq.size()-1] = deltaT[celli];
800 
801  tabulation_.add
802  (
803  phiq,
804  Rphiq,
805  mechRed_.nActiveSpecies(),
806  celli,
807  deltaT[celli]
808  );
809  }
810 
811  // When operations are done and if mechanism reduction is active,
812  // the number of species (which also affects nEqns) is set back
813  // to the total number of species (stored in the mechRed object)
814  if (reduction_)
815  {
816  setNSpecie(mechRed_.nSpecie());
817  }
818 
819  deltaTMin = min(deltaTChem_[celli], deltaTMin);
820  deltaTChem_[celli] = min(deltaTChem_[celli], deltaTChemMax_);
821  }
822 
823  // Set the RR vector (used in the solver)
824  for (label i=0; i<nSpecie_; i++)
825  {
826  RR_[i][celli] = rho0*(Y_[i] - Y0[i])/deltaT[celli];
827  }
828 
829  if (cpuLoad_)
830  {
831  chemistryCpuLoad.cpuTimeIncrement(celli);
832  }
833  }
834 
835  if (log_)
836  {
837  cpuSolveFile_()
838  << this->time().userTimeValue()
839  << " " << totalSolveCpuTime << endl;
840  }
841 
842  mechRed_.update();
843  tabulation_.update();
844 
845  if (reduction_)
846  {
847  this->thermo().syncSpeciesActive();
848  }
849 
850  return deltaTMin;
851 }
852 
853 
854 template<class ThermoType>
856 (
857  const scalar deltaT
858 )
859 {
860  // Don't allow the time-step to change more than a factor of 2
861  return min
862  (
864  2*deltaT
865  );
866 }
867 
868 
869 template<class ThermoType>
871 (
872  const scalarField& deltaT
873 )
874 {
875  return this->solve<scalarField>(deltaT);
876 }
877 
878 
879 template<class ThermoType>
882 {
884  (
886  (
887  "tc",
888  this->mesh(),
889  dimensionedScalar(dimTime, small),
890  extrapolatedCalculatedFvPatchScalarField::typeName
891  )
892  );
893  scalarField& tc = ttc.ref();
894 
895  if (!this->chemistry_)
896  {
897  ttc.ref().correctBoundaryConditions();
898  return ttc;
899  }
900 
901  tmp<volScalarField> trhovf(this->thermo().rho());
902  const volScalarField& rhovf = trhovf();
903 
904  const volScalarField& Tvf = this->thermo().T();
905  const volScalarField& pvf = this->thermo().p();
906 
907  reactionEvaluationScope scope(*this);
908 
909  forAll(rhovf, celli)
910  {
911  const scalar rho = rhovf[celli];
912  const scalar T = Tvf[celli];
913  const scalar p = pvf[celli];
914 
915  for (label i=0; i<nSpecie_; i++)
916  {
917  c_[i] = rho*Yvf_[i][celli]/specieThermos_[i].W();
918  }
919 
920  // A reaction's rate scale is calculated as its molar
921  // production rate divided by the total number of moles in the
922  // system.
923  //
924  // The system rate scale is the average of the reactions' rate
925  // scales weighted by the reactions' molar production rates. This
926  // weighting ensures that dominant reactions provide the largest
927  // contribution to the system rate scale.
928  //
929  // The system time scale is then the reciprocal of the system rate
930  // scale.
931  //
932  // Contributions from forward and reverse reaction rates are
933  // handled independently and identically so that reversible
934  // reactions produce the same result as the equivalent pair of
935  // irreversible reactions.
936 
937  scalar sumW = 0, sumWRateByCTot = 0;
938  forAll(reactions_, i)
939  {
940  const Reaction<ThermoType>& R = reactions_[i];
941  scalar omegaf, omegar;
942  R.omega(p, T, c_, celli, omegaf, omegar);
943 
944  scalar wf = 0;
945  forAll(R.rhs(), s)
946  {
947  wf += R.rhs()[s].stoichCoeff*omegaf;
948  }
949  sumW += wf;
950  sumWRateByCTot += sqr(wf);
951 
952  scalar wr = 0;
953  forAll(R.lhs(), s)
954  {
955  wr += R.lhs()[s].stoichCoeff*omegar;
956  }
957  sumW += wr;
958  sumWRateByCTot += sqr(wr);
959  }
960 
961  tc[celli] =
962  sumWRateByCTot == 0 ? vGreat : sumW/sumWRateByCTot*sum(c_);
963  }
964 
965  ttc.ref().correctBoundaryConditions();
966  return ttc;
967 }
968 
969 
970 template<class ThermoType>
973 {
974  tmp<volScalarField> tQdot
975  (
977  (
978  "Qdot",
979  this->mesh_,
981  )
982  );
983 
984  if (!this->chemistry_)
985  {
986  return tQdot;
987  }
988 
989  reactionEvaluationScope scope(*this);
990 
991  scalarField& Qdot = tQdot.ref();
992 
993  forAll(Yvf_, i)
994  {
995  forAll(Qdot, celli)
996  {
997  const scalar hi = specieThermos_[i].hf();
998  Qdot[celli] -= hi*RR_[i][celli];
999  }
1000  }
1001 
1002  return tQdot;
1003 }
1004 
1005 
1006 // ************************************************************************* //
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:433
Field with dimensions and associated with geometry type GeoMesh which is used to size the field and a...
static tmp< DimensionedField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Mesh &mesh, const dimensionSet &, const PrimitiveField< Type > &)
Return a temporary field constructed from name, mesh,.
Generic GeometricField class.
static tmp< GeometricField< Type, GeoMesh, PrimitiveField > > New(const word &name, const Internal &, const PtrList< Patch > &, 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:315
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
const Field0Type & oldTime() const
Return the old-time field.
Definition: OldTimeField.C:322
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 tmp< volScalarField::Internal > reactionRR(const label reactioni) const
Return the rate of reactioni [kmol/m^3/s].
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 > specieReactionRR(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 setSpecieInactive(const label speciei) const
Set specie inactive.
virtual void syncSpeciesActive() const =0
Synchronise the specie active flags.
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:197
Templated form of IOobject providing type information for file reading and header type checking.
Definition: IOobject.H:530
bool headerOk()
Read header (uses typeGlobalFile to find file) and check.
Foam::fvMesh mesh(Foam::IOobject(regionName, runTime.name(), runTime, Foam::IOobject::MUST_READ), false)
gmvFile<< "tracers "<< particles.size()<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().x()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().y()<< " ";}gmvFile<< nl;forAllConstIter(lagrangian::Cloud< passiveParticle >, particles, iter){ gmvFile<< iter().position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.name(), lagrangian::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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:258
void T(LagrangianPatchField< Type > &f, const LagrangianPatchField< Type > &f1)
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
const dimensionSet dimMoles
dimensioned< Type > sum(const DimensionedField< Type, GeoMesh, PrimitiveField > &df)
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
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