TDACChemistryModel.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-2019 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 "TDACChemistryModel.H"
27 #include "UniformField.H"
28 #include "localEulerDdtScheme.H"
29 #include "clockTime.H"
30 
31 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32 
33 template<class ReactionThermo, class ThermoType>
35 (
36  ReactionThermo& thermo
37 )
38 :
40  variableTimeStep_
41  (
42  this->mesh().time().controlDict().lookupOrDefault
43  (
44  "adjustTimeStep",
45  false
46  )
47  || fv::localEulerDdt::enabled(this->mesh())
48  ),
49  timeSteps_(0),
50  NsDAC_(this->nSpecie_),
51  completeC_(this->nSpecie_, 0),
52  reactionsDisabled_(this->reactions_.size(), false),
53  specieComp_(this->nSpecie_),
54  completeToSimplifiedIndex_(this->nSpecie_, -1),
55  simplifiedToCompleteIndex_(this->nSpecie_),
56  tabulationResults_
57  (
58  IOobject
59  (
60  thermo.phasePropertyName("TabulationResults"),
61  this->time().timeName(),
62  this->mesh(),
63  IOobject::NO_READ,
64  IOobject::AUTO_WRITE
65  ),
66  this->mesh(),
67  scalar(0)
68  )
69 {
70  basicSpecieMixture& composition = this->thermo().composition();
71 
72  // Store the species composition according to the species index
73  speciesTable speciesTab = composition.species();
74 
75  const HashTable<List<specieElement>>& specComp =
76  dynamicCast<const reactingMixture<ThermoType>&>(this->thermo())
77  .specieComposition();
78 
79  forAll(specieComp_, i)
80  {
81  specieComp_[i] = specComp[this->Y()[i].member()];
82  }
83 
85  (
86  *this,
87  *this
88  );
89 
90  // When the mechanism reduction method is used, the 'active' flag for every
91  // species should be initialized (by default 'active' is true)
92  if (mechRed_->active())
93  {
94  forAll(this->Y(), i)
95  {
96  IOobject header
97  (
98  this->Y()[i].name(),
99  this->mesh().time().timeName(),
100  this->mesh(),
101  IOobject::NO_READ
102  );
103 
104  // Check if the species file is provided, if not set inactive
105  // and NO_WRITE
106  if (!header.typeHeaderOk<volScalarField>(true))
107  {
108  composition.setInactive(i);
109  this->Y()[i].writeOpt() = IOobject::NO_WRITE;
110  }
111  }
112  }
113 
115  (
116  *this,
117  *this
118  );
119 
120  if (mechRed_->log())
121  {
122  cpuReduceFile_ = logFile("cpu_reduce.out");
123  nActiveSpeciesFile_ = logFile("nActiveSpecies.out");
124  }
125 
126  if (tabulation_->log())
127  {
128  cpuAddFile_ = logFile("cpu_add.out");
129  cpuGrowFile_ = logFile("cpu_grow.out");
130  cpuRetrieveFile_ = logFile("cpu_retrieve.out");
131  }
132 
133  if (mechRed_->log() || tabulation_->log())
134  {
135  cpuSolveFile_ = logFile("cpu_solve.out");
136  }
137 }
138 
139 
140 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
141 
142 template<class ReactionThermo, class ThermoType>
144 {}
145 
146 
147 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
148 
149 template<class ReactionThermo, class ThermoType>
151 (
152  const scalarField& c, // Contains all species even when mechRed is active
153  const scalar T,
154  const scalar p,
155  scalarField& dcdt
156 ) const
157 {
158  const bool reduced = mechRed_->active();
159 
160  scalar pf, cf, pr, cr;
161  label lRef, rRef;
162 
163  dcdt = Zero;
164 
165  forAll(this->reactions_, i)
166  {
167  if (!reactionsDisabled_[i])
168  {
169  const Reaction<ThermoType>& R = this->reactions_[i];
170 
171  scalar omegai = R.omega
172  (
173  p, T, c, pf, cf, lRef, pr, cr, rRef
174  );
175 
176  forAll(R.lhs(), s)
177  {
178  label si = R.lhs()[s].index;
179  if (reduced)
180  {
181  si = completeToSimplifiedIndex_[si];
182  }
183 
184  const scalar sl = R.lhs()[s].stoichCoeff;
185  dcdt[si] -= sl*omegai;
186  }
187  forAll(R.rhs(), s)
188  {
189  label si = R.rhs()[s].index;
190  if (reduced)
191  {
192  si = completeToSimplifiedIndex_[si];
193  }
194 
195  const scalar sr = R.rhs()[s].stoichCoeff;
196  dcdt[si] += sr*omegai;
197  }
198  }
199  }
200 }
201 
202 
203 template<class ReactionThermo, class ThermoType>
205 (
206  const Reaction<ThermoType>& R,
207  const scalarField& c, // Contains all species even when mechRed is active
208  const scalar T,
209  const scalar p,
210  scalar& pf,
211  scalar& cf,
212  label& lRef,
213  scalar& pr,
214  scalar& cr,
215  label& rRef
216 ) const
217 {
218  const scalar kf = R.kf(p, T, c);
219  const scalar kr = R.kr(kf, p, T, c);
220 
221  const label Nl = R.lhs().size();
222  const label Nr = R.rhs().size();
223 
224  label slRef = 0;
225  lRef = R.lhs()[slRef].index;
226 
227  pf = kf;
228  for (label s=1; s<Nl; s++)
229  {
230  const label si = R.lhs()[s].index;
231 
232  if (c[si] < c[lRef])
233  {
234  const scalar exp = R.lhs()[slRef].exponent;
235  pf *= pow(max(c[lRef], 0), exp);
236  lRef = si;
237  slRef = s;
238  }
239  else
240  {
241  const scalar exp = R.lhs()[s].exponent;
242  pf *= pow(max(c[si], 0), exp);
243  }
244  }
245  cf = max(c[lRef], 0);
246 
247  {
248  const scalar exp = R.lhs()[slRef].exponent;
249  if (exp < 1)
250  {
251  if (cf > small)
252  {
253  pf *= pow(cf, exp - 1);
254  }
255  else
256  {
257  pf = 0;
258  }
259  }
260  else
261  {
262  pf *= pow(cf, exp - 1);
263  }
264  }
265 
266  label srRef = 0;
267  rRef = R.rhs()[srRef].index;
268 
269  // Find the matrix element and element position for the rhs
270  pr = kr;
271  for (label s=1; s<Nr; s++)
272  {
273  const label si = R.rhs()[s].index;
274  if (c[si] < c[rRef])
275  {
276  const scalar exp = R.rhs()[srRef].exponent;
277  pr *= pow(max(c[rRef], 0), exp);
278  rRef = si;
279  srRef = s;
280  }
281  else
282  {
283  const scalar exp = R.rhs()[s].exponent;
284  pr *= pow(max(c[si], 0), exp);
285  }
286  }
287  cr = max(c[rRef], 0);
288 
289  {
290  const scalar exp = R.rhs()[srRef].exponent;
291  if (exp < 1)
292  {
293  if (cr > small)
294  {
295  pr *= pow(cr, exp - 1);
296  }
297  else
298  {
299  pr = 0;
300  }
301  }
302  else
303  {
304  pr *= pow(cr, exp - 1);
305  }
306  }
307 
308  return pf*cf - pr*cr;
309 }
310 
311 
312 template<class ReactionThermo, class ThermoType>
314 (
315  const scalar time,
316  const scalarField& c,
317  scalarField& dcdt
318 ) const
319 {
320  const bool reduced = mechRed_->active();
321 
322  const scalar T = c[this->nSpecie_];
323  const scalar p = c[this->nSpecie_ + 1];
324 
325  if (reduced)
326  {
327  // When using DAC, the ODE solver submit a reduced set of species the
328  // complete set is used and only the species in the simplified mechanism
329  // are updated
330  this->c_ = completeC_;
331 
332  // Update the concentration of the species in the simplified mechanism
333  // the other species remain the same and are used only for third-body
334  // efficiencies
335  for (label i=0; i<NsDAC_; i++)
336  {
337  this->c_[simplifiedToCompleteIndex_[i]] = max(c[i], 0);
338  }
339  }
340  else
341  {
342  for (label i=0; i<this->nSpecie(); i++)
343  {
344  this->c_[i] = max(c[i], 0);
345  }
346  }
347 
348  omega(this->c_, T, p, dcdt);
349 
350  // Constant pressure
351  // dT/dt = ...
352  scalar rho = 0;
353  for (label i=0; i<this->c_.size(); i++)
354  {
355  const scalar W = this->specieThermo_[i].W();
356  rho += W*this->c_[i];
357  }
358 
359  scalar cp = 0;
360  for (label i=0; i<this->c_.size(); i++)
361  {
362  // cp function returns [J/kmol/K]
363  cp += this->c_[i]*this->specieThermo_[i].cp(p, T);
364  }
365  cp /= rho;
366 
367  // When mechanism reduction is active
368  // dT is computed on the reduced set since dcdt is null
369  // for species not involved in the simplified mechanism
370  scalar dT = 0;
371  for (label i=0; i<this->nSpecie_; i++)
372  {
373  label si;
374  if (reduced)
375  {
376  si = simplifiedToCompleteIndex_[i];
377  }
378  else
379  {
380  si = i;
381  }
382 
383  // ha function returns [J/kmol]
384  const scalar hi = this->specieThermo_[si].ha(p, T);
385  dT += hi*dcdt[i];
386  }
387  dT /= rho*cp;
388 
389  dcdt[this->nSpecie_] = -dT;
390 
391  // dp/dt = ...
392  dcdt[this->nSpecie_ + 1] = 0;
393 }
394 
395 
396 template<class ReactionThermo, class ThermoType>
398 (
399  const scalar t,
400  const scalarField& c,
401  scalarField& dcdt,
403 ) const
404 {
405  const bool reduced = mechRed_->active();
406 
407  // If the mechanism reduction is active, the computed Jacobian
408  // is compact (size of the reduced set of species)
409  // but according to the information of the complete set
410  // (i.e. for the third-body efficiencies)
411 
412  const scalar T = c[this->nSpecie_];
413  const scalar p = c[this->nSpecie_ + 1];
414 
415  if (reduced)
416  {
417  this->c_ = completeC_;
418  for (label i=0; i<NsDAC_; i++)
419  {
420  this->c_[simplifiedToCompleteIndex_[i]] = max(c[i], 0);
421  }
422  }
423  else
424  {
425  forAll(this->c_, i)
426  {
427  this->c_[i] = max(c[i], 0);
428  }
429  }
430 
431  J = Zero;
432  dcdt = Zero;
433  scalarField hi(this->c_.size());
434  scalarField cpi(this->c_.size());
435  forAll(hi, i)
436  {
437  hi[i] = this->specieThermo_[i].ha(p, T);
438  cpi[i] = this->specieThermo_[i].cp(p, T);
439  }
440 
441  scalar omegaI = 0;
442 
443  forAll(this->reactions_, ri)
444  {
445  if (!reactionsDisabled_[ri])
446  {
447  const Reaction<ThermoType>& R = this->reactions_[ri];
448  scalar kfwd, kbwd;
449  R.dwdc
450  (
451  p,
452  T,
453  this->c_,
454  J,
455  dcdt,
456  omegaI,
457  kfwd,
458  kbwd,
459  reduced,
460  completeToSimplifiedIndex_
461  );
462  R.dwdT
463  (
464  p,
465  T,
466  this->c_,
467  omegaI,
468  kfwd,
469  kbwd,
470  J,
471  reduced,
472  completeToSimplifiedIndex_,
473  this->nSpecie_
474  );
475  }
476  }
477 
478  // The species derivatives of the temperature term are partially computed
479  // while computing dwdc, they are completed hereunder:
480  scalar cpMean = 0;
481  scalar dcpdTMean = 0;
482  forAll(this->c_, i)
483  {
484  cpMean += this->c_[i]*cpi[i]; // J/(m^3 K)
485  // Already multiplied by rho
486  dcpdTMean += this->c_[i]*this->specieThermo_[i].dcpdT(p, T);
487  }
488 
489  scalar dTdt = 0;
490  forAll(hi, i)
491  {
492  if (reduced)
493  {
494  const label si = completeToSimplifiedIndex_[i];
495  if (si != -1)
496  {
497  dTdt += hi[i]*dcdt[si]; // J/(m^3 s)
498  }
499  }
500  else
501  {
502  dTdt += hi[i]*dcdt[i]; // J/(m^3 s)
503  }
504  }
505  dTdt /= -cpMean; // K/s
506  dcdt[this->nSpecie_] = dTdt;
507 
508  for (label i = 0; i < this->nSpecie_; i++)
509  {
510  J(this->nSpecie_, i) = 0;
511  for (label j = 0; j < this->nSpecie_; j++)
512  {
513  const label sj = reduced ? simplifiedToCompleteIndex_[j] : j;
514  J(this->nSpecie_, i) += hi[sj]*J(j, i);
515  }
516  const label si = reduced ? simplifiedToCompleteIndex_[i] : i;
517  J(this->nSpecie_, i) += cpi[si]*dTdt; // J/(mol s)
518  J(this->nSpecie_, i) /= -cpMean; // K/s / (mol/m^3)
519  }
520 
521  // ddT of dTdt
522  J(this->nSpecie_, this->nSpecie_) = 0;
523  for (label i = 0; i < this->nSpecie_; i++)
524  {
525  const label si = reduced ? simplifiedToCompleteIndex_[i] : i;
526  J(this->nSpecie_, this->nSpecie_) +=
527  cpi[si]*dcdt[i]
528  + hi[si]*J(i, this->nSpecie_);
529  }
530  J(this->nSpecie_, this->nSpecie_) += dTdt*dcpdTMean;
531  J(this->nSpecie_, this->nSpecie_) /= -cpMean;
532  J(this->nSpecie_, this->nSpecie_) += dTdt/T;
533 }
534 
535 
536 template<class ReactionThermo, class ThermoType>
537 template<class DeltaTType>
539 (
540  const DeltaTType& deltaT
541 )
542 {
543  // Increment counter of time-step
544  timeSteps_++;
545 
546  const bool reduced = mechRed_->active();
547 
548  label nAdditionalEqn = (tabulation_->variableTimeStep() ? 1 : 0);
549 
550  basicSpecieMixture& composition = this->thermo().composition();
551 
552  // CPU time analysis
553  const clockTime clockTime_= clockTime();
554  clockTime_.timeIncrement();
555  scalar reduceMechCpuTime_ = 0;
556  scalar addNewLeafCpuTime_ = 0;
557  scalar growCpuTime_ = 0;
558  scalar solveChemistryCpuTime_ = 0;
559  scalar searchISATCpuTime_ = 0;
560 
561  this->resetTabulationResults();
562 
563  // Average number of active species
564  scalar nActiveSpecies = 0;
565  scalar nAvg = 0;
566 
568 
569  scalar deltaTMin = great;
570 
571  if (!this->chemistry_)
572  {
573  return deltaTMin;
574  }
575 
576  const volScalarField rho
577  (
578  IOobject
579  (
580  "rho",
581  this->time().timeName(),
582  this->mesh(),
583  IOobject::NO_READ,
584  IOobject::NO_WRITE,
585  false
586  ),
587  this->thermo().rho()
588  );
589 
590  const scalarField& T = this->thermo().T();
591  const scalarField& p = this->thermo().p();
592 
593  scalarField c(this->nSpecie_);
594  scalarField c0(this->nSpecie_);
595 
596  // Composition vector (Yi, T, p)
597  scalarField phiq(this->nEqns() + nAdditionalEqn);
598 
599  scalarField Rphiq(this->nEqns() + nAdditionalEqn);
600 
601  forAll(rho, celli)
602  {
603  const scalar rhoi = rho[celli];
604  scalar pi = p[celli];
605  scalar Ti = T[celli];
606 
607  for (label i=0; i<this->nSpecie_; i++)
608  {
609  c[i] = rhoi*this->Y_[i][celli]/this->specieThermo_[i].W();
610  c0[i] = c[i];
611  phiq[i] = this->Y()[i][celli];
612  }
613  phiq[this->nSpecie()]=Ti;
614  phiq[this->nSpecie() + 1]=pi;
615  if (tabulation_->variableTimeStep())
616  {
617  phiq[this->nSpecie() + 2] = deltaT[celli];
618  }
619 
620 
621  // Initialise time progress
622  scalar timeLeft = deltaT[celli];
623 
624  // Not sure if this is necessary
625  Rphiq = Zero;
626 
627  clockTime_.timeIncrement();
628 
629  // When tabulation is active (short-circuit evaluation for retrieve)
630  // It first tries to retrieve the solution of the system with the
631  // information stored through the tabulation method
632  if (tabulation_->active() && tabulation_->retrieve(phiq, Rphiq))
633  {
634  // Retrieved solution stored in Rphiq
635  for (label i=0; i<this->nSpecie(); i++)
636  {
637  c[i] = rhoi*Rphiq[i]/this->specieThermo_[i].W();
638  }
639 
640  searchISATCpuTime_ += clockTime_.timeIncrement();
641  }
642  // This position is reached when tabulation is not used OR
643  // if the solution is not retrieved.
644  // In the latter case, it adds the information to the tabulation
645  // (it will either expand the current data or add a new stored point).
646  else
647  {
648  // Store total time waiting to attribute to add or grow
649  scalar timeTmp = clockTime_.timeIncrement();
650 
651  if (reduced)
652  {
653  // Reduce mechanism change the number of species (only active)
654  mechRed_->reduceMechanism(c, Ti, pi);
655  nActiveSpecies += mechRed_->NsSimp();
656  nAvg++;
657  scalar timeIncr = clockTime_.timeIncrement();
658  reduceMechCpuTime_ += timeIncr;
659  timeTmp += timeIncr;
660  }
661 
662  // Calculate the chemical source terms
663  while (timeLeft > small)
664  {
665  scalar dt = timeLeft;
666  if (reduced)
667  {
668  // completeC_ used in the overridden ODE methods
669  // to update only the active species
670  completeC_ = c;
671 
672  // Solve the reduced set of ODE
673  this->solve
674  (
675  simplifiedC_, Ti, pi, dt, this->deltaTChem_[celli]
676  );
677 
678  for (label i=0; i<NsDAC_; i++)
679  {
680  c[simplifiedToCompleteIndex_[i]] = simplifiedC_[i];
681  }
682  }
683  else
684  {
685  this->solve(c, Ti, pi, dt, this->deltaTChem_[celli]);
686  }
687  timeLeft -= dt;
688  }
689 
690  {
691  scalar timeIncr = clockTime_.timeIncrement();
692  solveChemistryCpuTime_ += timeIncr;
693  timeTmp += timeIncr;
694  }
695 
696  // If tabulation is used, we add the information computed here to
697  // the stored points (either expand or add)
698  if (tabulation_->active())
699  {
700  forAll(c, i)
701  {
702  Rphiq[i] = c[i]/rhoi*this->specieThermo_[i].W();
703  }
704  if (tabulation_->variableTimeStep())
705  {
706  Rphiq[Rphiq.size()-3] = Ti;
707  Rphiq[Rphiq.size()-2] = pi;
708  Rphiq[Rphiq.size()-1] = deltaT[celli];
709  }
710  else
711  {
712  Rphiq[Rphiq.size()-2] = Ti;
713  Rphiq[Rphiq.size()-1] = pi;
714  }
715  label growOrAdd =
716  tabulation_->add(phiq, Rphiq, rhoi, deltaT[celli]);
717  if (growOrAdd)
718  {
719  this->setTabulationResultsAdd(celli);
720  addNewLeafCpuTime_ += clockTime_.timeIncrement() + timeTmp;
721  }
722  else
723  {
724  this->setTabulationResultsGrow(celli);
725  growCpuTime_ += clockTime_.timeIncrement() + timeTmp;
726  }
727  }
728 
729  // When operations are done and if mechanism reduction is active,
730  // the number of species (which also affects nEqns) is set back
731  // to the total number of species (stored in the mechRed object)
732  if (reduced)
733  {
734  this->nSpecie_ = mechRed_->nSpecie();
735  }
736  deltaTMin = min(this->deltaTChem_[celli], deltaTMin);
737 
738  this->deltaTChem_[celli] =
739  min(this->deltaTChem_[celli], this->deltaTChemMax_);
740  }
741 
742  // Set the RR vector (used in the solver)
743  for (label i=0; i<this->nSpecie_; i++)
744  {
745  this->RR_[i][celli] =
746  (c[i] - c0[i])*this->specieThermo_[i].W()/deltaT[celli];
747  }
748  }
749 
750  if (mechRed_->log() || tabulation_->log())
751  {
752  cpuSolveFile_()
753  << this->time().timeOutputValue()
754  << " " << solveChemistryCpuTime_ << endl;
755  }
756 
757  if (mechRed_->log())
758  {
759  cpuReduceFile_()
760  << this->time().timeOutputValue()
761  << " " << reduceMechCpuTime_ << endl;
762  }
763 
764  if (tabulation_->active())
765  {
766  // Every time-step, look if the tabulation should be updated
767  tabulation_->update();
768 
769  // Write the performance of the tabulation
770  tabulation_->writePerformance();
771 
772  if (tabulation_->log())
773  {
774  cpuRetrieveFile_()
775  << this->time().timeOutputValue()
776  << " " << searchISATCpuTime_ << endl;
777 
778  cpuGrowFile_()
779  << this->time().timeOutputValue()
780  << " " << growCpuTime_ << endl;
781 
782  cpuAddFile_()
783  << this->time().timeOutputValue()
784  << " " << addNewLeafCpuTime_ << endl;
785  }
786  }
787 
788  if (reduced && nAvg && mechRed_->log())
789  {
790  // Write average number of species
791  nActiveSpeciesFile_()
792  << this->time().timeOutputValue()
793  << " " << nActiveSpecies/nAvg << endl;
794  }
795 
796  if (Pstream::parRun())
797  {
798  List<bool> active(composition.active());
799  Pstream::listCombineGather(active, orEqOp<bool>());
800  Pstream::listCombineScatter(active);
801 
802  forAll(active, i)
803  {
804  if (active[i])
805  {
806  composition.setActive(i);
807  }
808  }
809  }
810 
811  forAll(this->Y(), i)
812  {
813  if (composition.active(i))
814  {
815  this->Y()[i].writeOpt() = IOobject::AUTO_WRITE;
816  }
817  }
818 
819  return deltaTMin;
820 }
821 
822 
823 template<class ReactionThermo, class ThermoType>
825 (
826  const scalar deltaT
827 )
828 {
829  // Don't allow the time-step to change more than a factor of 2
830  return min
831  (
833  2*deltaT
834  );
835 }
836 
837 
838 template<class ReactionThermo, class ThermoType>
840 (
841  const scalarField& deltaT
842 )
843 {
844  return this->solve<scalarField>(deltaT);
845 }
846 
847 
848 template<class ReactionThermo, class ThermoType>
851 (
852  const label celli
853 )
854 {
855  tabulationResults_[celli] = 0;
856 }
857 
858 
859 template<class ReactionThermo, class ThermoType>
862 {
863  tabulationResults_[celli] = 1;
864 }
865 
866 
867 template<class ReactionThermo, class ThermoType>
870 (
871  const label celli
872 )
873 {
874  tabulationResults_[celli] = 2;
875 }
876 
877 
878 // ************************************************************************* //
A class representing the concept of a uniform field which stores only the single value and providing ...
Definition: UniformField.H:47
Extends StandardChemistryModel by adding the TDAC method.
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
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
bool typeHeaderOk(const bool checkType=true)
Read header (uses typeFilePath to find file) and check header.
bool active(label speciei) const
Return true for active species.
void setTabulationResultsGrow(const label celli)
void setTabulationResultsRetrieve(const label celli)
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void dwdT(const scalar p, const scalar T, const scalarField &c, const scalar omegaI, const scalar kfwd, const scalar kbwd, scalarSquareMatrix &J, const bool reduced, const List< label > &c2s, const label indexT) const
Derivative of the net reaction rate for each species involved.
Definition: Reaction.C:537
basicSpecieMixture & composition
TDACChemistryModel(ReactionThermo &thermo)
Construct from thermo.
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:163
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:256
virtual ~TDACChemistryModel()
Destructor.
const speciesTable & species() const
Return the table of species.
const List< specieCoeffs > & lhs() const
Return the components of the left hand side.
Definition: ReactionI.H:57
virtual scalar kr(const scalar kfwd, const scalar p, const scalar T, const scalarField &c) const =0
Reverse rate constant from the given forward rate constant.
rhoReactionThermo & thermo
Definition: createFields.H:28
const label nSpecie
void omega(const scalar p, const scalar T, const scalarField &c, scalarField &dcdt) const
Net reaction rate for individual species.
Definition: Reaction.C:228
Specialization of basicMultiComponentMixture for a mixture consisting of a number for molecular speci...
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
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
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:56
tmp< DimensionedField< TypeR, GeoMesh > > New(const tmp< DimensionedField< TypeR, GeoMesh >> &tdf1, const word &name, const dimensionSet &dimensions)
void setInactive(label speciei)
Set speciei inactive.
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.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
dynamicFvMesh & mesh
dimensionedScalar exp(const dimensionedScalar &ds)
virtual void jacobian(const scalar t, const scalarField &c, scalarField &dcdt, scalarSquareMatrix &J) const
Calculate the Jacobian of the system.
rhoEqn solve()
word timeName
Definition: getTimeIndex.H:3
static const zero Zero
Definition: zero.H:97
void dwdc(const scalar p, const scalar T, const scalarField &c, scalarSquareMatrix &J, scalarField &dcdt, scalar &omegaI, scalar &kfwd, scalar &kbwd, const bool reduced, const List< label > &c2s) const
Derivative of the net reaction rate for each species involved.
Definition: Reaction.C:374
An STL-conforming hash table.
Definition: HashTable.H:61
void setTabulationResultsAdd(const label celli)
void setActive(label speciei)
Set speciei active.
const volScalarField & cp
virtual void derivatives(const scalar t, const scalarField &c, scalarField &dcdt) const
Calculate the derivatives in dydx.
Info<< "Predicted p max-min : "<< max(p).value()<< " "<< min(p).value()<< endl;rho==max(rho0+psi *p, rhoMin);# 1 "/home/ubuntu/OpenFOAM-7/applications/solvers/multiphase/cavitatingFoam/alphavPsi.H" 1{ alphav=max(min((rho - rholSat)/(rhovSat - rholSat), scalar(1)), scalar(0));alphal=1.0 - alphav;Info<< "max-min alphav: "<< max(alphav).value()<< " "<< min(alphav).value()<< endl;psiModel-> correct()
Definition: pEqn.H:68
dimensioned< Type > min(const dimensioned< Type > &, const dimensioned< Type > &)
const volScalarField & T
word name(const complex &)
Return a string representation of a complex.
Definition: complex.C:47
dimensionedScalar pow(const dimensionedScalar &ds, const dimensionedScalar &expt)
bool active(const label i) const
A wordList with hashed indices for faster lookup by name.
const scalarList W(::W(thermo))
PtrList< volScalarField > & Y
const dimensionedScalar c
Speed of light in a vacuum.
double timeIncrement() const
Return time (in seconds) since last call to timeIncrement()
Definition: clockTime.C:62
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
Definition: ReactionI.H:64
IOobject defines the attributes of an object for which implicit objectRegistry management is supporte...
Definition: IOobject.H:92
virtual scalar kf(const scalar p, const scalar T, const scalarField &c) const =0
Forward rate constant.
Starts timing (using rtc) and returns elapsed time from start. Better resolution (2uSec instead of ~2...
Definition: clockTime.H:49