DAC.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-2021 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 "DAC.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ThermoType>
32 (
33  const IOdictionary& dict,
35 )
36 :
38  searchInitSet_(),
39  zprime_(0),
40  nbCLarge_(3),
41  sC_(this->nSpecie_,0),
42  sH_(this->nSpecie_,0),
43  sO_(this->nSpecie_,0),
44  sN_(this->nSpecie_,0),
45  CO2Id_(-1),
46  COId_(-1),
47  HO2Id_(-1),
48  H2OId_(-1),
49  NOId_(-1),
50  automaticSIS_(true),
51  phiTol_(this->tolerance()),
52  NOxThreshold_(1800),
53  CO2Name_
54  (
55  dict.subDict("reduction").lookupOrDefault<word>
56  (
57  "CO2Name","CO2"
58  )
59  ),
60  COName_
61  (
62  dict.subDict("reduction").lookupOrDefault<word>
63  (
64  "COName","CO"
65  )
66  ),
67  HO2Name_
68  (
69  dict.subDict("reduction").lookupOrDefault<word>
70  (
71  "HO2Name","HO2"
72  )
73  ),
74  H2OName_
75  (
76  dict.subDict("reduction").lookupOrDefault<word>
77  (
78  "H2OName","H2O"
79  )
80  ),
81  NOName_
82  (
83  dict.subDict("reduction").lookupOrDefault<word>
84  (
85  "NOName","NO"
86  )
87  ),
88  forceFuelInclusion_(false)
89 {
90  const wordHashSet initSet(this->coeffsDict_.lookup("initialSet"));
91  forAllConstIter(wordHashSet, initSet, iter)
92  {
93  searchInitSet_.append(chemistry.mixture().species()[iter.key()]);
94  }
95 
96  if (this->coeffsDict_.found("automaticSIS"))
97  {
98  automaticSIS_.readIfPresent("automaticSIS", this->coeffsDict_);
99  }
100 
101  if (this->coeffsDict_.found("forceFuelInclusion"))
102  {
103  forceFuelInclusion_.readIfPresent
104  (
105  "forceFuelInclusion",this->coeffsDict_
106  );
107  }
108 
109  if (this->coeffsDict_.found("phiTol"))
110  {
111  phiTol_ = this->coeffsDict_.template lookup<scalar>("phiTol");
112  }
113 
114  if (this->coeffsDict_.found("NOxThreshold"))
115  {
116  NOxThreshold_ =
117  this->coeffsDict_.template lookup<scalar>("NOxThreshold");
118  }
119 
120  for (label i=0; i<this->nSpecie_; i++)
121  {
122  const List<specieElement>& curSpecieComposition =
123  chemistry.mixture().specieComposition(i);
124 
125  // For all elements in the current species
126  forAll(curSpecieComposition, j)
127  {
128  const specieElement& curElement =
129  curSpecieComposition[j];
130  if (curElement.name() == "C")
131  {
132  sC_[i] = curElement.nAtoms();
133  }
134  else if (curElement.name() == "H")
135  {
136  sH_[i] = curElement.nAtoms();
137  }
138  else if (curElement.name() == "O")
139  {
140  sO_[i] = curElement.nAtoms();
141  }
142  else if (curElement.name() == "N")
143  {
144  sN_[i] = curElement.nAtoms();
145  }
146  else
147  {
148  Info<< "element not considered"<<endl;
149  }
150  }
151  if (this->chemistry_.Y()[i].member() == CO2Name_)
152  {
153  CO2Id_ = i;
154  }
155  else if (this->chemistry_.Y()[i].member() == COName_)
156  {
157  COId_ = i;
158  }
159  else if (this->chemistry_.Y()[i].member() == HO2Name_)
160  {
161  HO2Id_ = i;
162  }
163  else if (this->chemistry_.Y()[i].member() == H2OName_)
164  {
165  H2OId_ = i;
166  }
167  else if (this->chemistry_.Y()[i].member() == NOName_)
168  {
169  NOId_ = i;
170  }
171  }
172 
173  if ((CO2Id_==-1 || COId_==-1 || HO2Id_==-1 || H2OId_==-1) && automaticSIS_)
174  {
176  << "The name of the species used in automatic SIS are not found in "
177  << " the mechanism. You should either set the name for CO2, CO, H2O"
178  << " and HO2 properly or set automaticSIS to off "
179  << exit(FatalError);
180  }
181 
182  // To compute zprime, the fuel species should be specified.
183  // According to the given mass fraction, an equivalent O/C ratio is computed
184  if (automaticSIS_)
185  {
186  List<Tuple2<word, scalar>> fuelSpeciesEntry
187  (
188  this->coeffsDict_.lookup("fuelSpecies")
189  );
190 
191  fuelSpecies_.setSize(fuelSpeciesEntry.size());
192  fuelSpeciesID_.setSize(fuelSpeciesEntry.size());
193  fuelSpeciesProp_.setSize(fuelSpeciesEntry.size());
194  scalar Mmtot(0.0);
195 
196  forAll(fuelSpeciesEntry, i)
197  {
198  fuelSpecies_[i] = fuelSpeciesEntry[i].first();
199  fuelSpeciesProp_[i] = fuelSpeciesEntry[i].second();
200  fuelSpeciesID_[i] =
201  this->chemistry_.mixture().species()[fuelSpecies_[i]];
202  scalar curMm =
203  this->chemistry_.specieThermos()[fuelSpeciesID_[i]].W();
204  Mmtot += fuelSpeciesProp_[i]/curMm;
205  }
206 
207  this->coeffsDict_.readIfPresent("nbCLarge", nbCLarge_);
208 
209  Mmtot = 1.0/Mmtot;
210  scalar nbC(0.0);
211  scalar nbO(0.0);
212  forAll(fuelSpecies_, i)
213  {
214  label curID = fuelSpeciesID_[i];
215  scalar curMm = this->chemistry_.specieThermos()[curID].W();
216 
217  nbC += fuelSpeciesProp_[i]*Mmtot/curMm*sC_[curID];
218  nbO += fuelSpeciesProp_[i]*Mmtot/curMm*sO_[curID];
219  }
220  zprime_ = nbO/nbC;
221  }
222 }
223 
224 
225 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
226 
227 template<class ThermoType>
229 {}
230 
231 
232 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
233 
234 template<class ThermoType>
236 (
237  const scalar p,
238  const scalar T,
239  const scalarField& c,
240  const label li
241 )
242 {
243  scalarField& completeC(this->chemistry_.completeC());
244  scalarField c1(this->chemistry_.nEqns(), 0.0);
245  for(label i=0; i<this->nSpecie_; i++)
246  {
247  c1[i] = c[i];
248  completeC[i] = c[i];
249  }
250 
251  c1[this->nSpecie_] = T;
252  c1[this->nSpecie_+1] = p;
253 
254  // Compute the rAB matrix
255  RectangularMatrix<scalar> rABNum(this->nSpecie_,this->nSpecie_,0.0);
256  scalarField PA(this->nSpecie_,0.0);
257  scalarField CA(this->nSpecie_,0.0);
258 
259  // Number of initialised rAB for each lines
260  Field<label> NbrABInit(this->nSpecie_,0);
261  // Position of the initialised rAB, -1 when not initialised
262  RectangularMatrix<label> rABPos(this->nSpecie_, this->nSpecie_, -1);
263  // Index of the other species involved in the rABNum
264  RectangularMatrix<label> rABOtherSpec(this->nSpecie_, this->nSpecie_, -1);
265 
266  forAll(this->chemistry_.reactions(), i)
267  {
268  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
269 
270  // for each reaction compute omegai
271  scalar omegaf, omegar;
272  const scalar omegai = R.omega(p, T, c1, li, omegaf, omegar);
273 
274  // Then for each pair of species composing this reaction,
275  // compute the rAB matrix (separate the numerator and
276  // denominator)
277 
278  // While computing the rAB for all the species involved in the reaction
279  // we should consider that one can write a reaction A+B=2C as A+B=C+C
280  // In this case, the following algorithm only take once the effect
281  // of the species. It stores the species encountered in the reaction but
282  // use another list to see if this species has already been used
283 
284  DynamicList<scalar> wA(R.lhs().size()+R.rhs().size());
285  DynamicList<label> wAID(R.lhs().size()+R.rhs().size());
286 
287  forAll(R.lhs(), s) // Compute rAB for all species in the left hand side
288  {
289  label ss = R.lhs()[s].index;
290  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
291  List<bool> deltaBi(this->nSpecie_, false);
292  FIFOStack<label> usedIndex;
293  forAll(R.lhs(), j)
294  {
295  label sj = R.lhs()[j].index;
296  usedIndex.push(sj);
297  deltaBi[sj] = true;
298  }
299  forAll(R.rhs(), j)
300  {
301  label sj = R.rhs()[j].index;
302  usedIndex.push(sj);
303  deltaBi[sj] = true;
304  }
305 
306  // Disable for self reference (by definition rAA=0)
307  deltaBi[ss] = false;
308 
309  while(!usedIndex.empty())
310  {
311  label curIndex = usedIndex.pop();
312  if (deltaBi[curIndex])
313  {
314  // Disable to avoid counting it more than once
315  deltaBi[curIndex] = false;
316  // Test if this rAB is not initialised
317  if (rABPos(ss, curIndex)==-1)
318  {
319  // It starts at rABPos(ss, sj)=0
320  rABPos(ss, curIndex) = NbrABInit[ss];
321  NbrABInit[ss]++;
322  // to avoid overflow
323  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
324  // store the other specie involved
325  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
326  }
327  else
328  {
329  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
330  }
331  }
332  }
333 
334  bool found(false);
335  forAll(wAID, id)
336  {
337  if (ss==wAID[id])
338  {
339  wA[id] += sl*omegai;
340  found = true;
341  }
342  }
343  if (!found)
344  {
345  wA.append(sl*omegai);
346  wAID.append(ss);
347  }
348  }
349 
350  forAll(R.rhs(), s) // Compute rAB for all species in the right hand side
351  {
352  label ss = R.rhs()[s].index;
353  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
354  List<bool> deltaBi(this->nSpecie_, false);
355  FIFOStack<label> usedIndex;
356  forAll(R.lhs(), j)
357  {
358  label sj = R.lhs()[j].index;
359  usedIndex.push(sj);
360  deltaBi[sj] = true;
361  }
362  forAll(R.rhs(), j)
363  {
364  label sj = R.rhs()[j].index;
365  usedIndex.push(sj);
366  deltaBi[sj] = true;
367  }
368 
369  // Disable for self reference (by definition rAA=0)
370  deltaBi[ss] = false;
371 
372  while(!usedIndex.empty())
373  {
374  label curIndex = usedIndex.pop();
375  if (deltaBi[curIndex])
376  {
377  // Disable to avoid counting it more than once
378  deltaBi[curIndex] = false;
379 
380  // Test if this rAB is not initialised
381  if (rABPos(ss, curIndex) == -1)
382  {
383  // it starts at rABPos(ss, sj)=0
384  rABPos(ss, curIndex) = NbrABInit[ss];
385  NbrABInit[ss]++;
386  rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
387  rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
388  }
389  else
390  {
391  rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
392  }
393  }
394  }
395 
396  bool found(false);
397  forAll(wAID, id)
398  {
399  if (ss==wAID[id])
400  {
401  wA[id] += sl*omegai;
402  found = true;
403  }
404  }
405  if (!found)
406  {
407  wA.append(sl*omegai);
408  wAID.append(ss);
409  }
410  }
411  wAID.shrink();
412 
413  // Now that every species of the reactions has been visited, we can
414  // compute the production and consumption rate. This way, it avoids
415  // getting wrong results when species are present in both lhs and rhs
416  forAll(wAID, id)
417  {
418  if (wA[id] > 0.0)
419  {
420  if (PA[wAID[id]] == 0.0)
421  {
422  PA[wAID[id]] = wA[id];
423  }
424  else
425  {
426  PA[wAID[id]] += wA[id];
427  }
428  }
429  else
430  {
431  if (CA[wAID[id]] == 0.0)
432  {
433  CA[wAID[id]] = -wA[id];
434  }
435  else
436  {
437  CA[wAID[id]] += -wA[id];
438  }
439  }
440  }
441  }
442  // rii = 0.0 by definition
443 
444  scalar phiLarge(0.0);
445  scalar phiProgress(0.0);
446  if (automaticSIS_)
447  {
448  // Compute the progress equivalence ratio
449  // and the equivalence ratio for fuel decomposition
450  label nElements = 4; // 4 main elements (C, H, O, N)
451 
452  // Total number of C, H and O (in this order)
453  scalarList Na(nElements,0.0);
454  scalarList Nal(nElements,0.0); // for large hydrocarbons
455 
456  for (label i=0; i<this->nSpecie_; i++)
457  {
458  // Complete combustion products are not considered
459  if
460  (
461  this->chemistry_.Y()[i].member() == "CO2"
462  || this->chemistry_.Y()[i].member() == "H2O"
463  )
464  {
465  continue;
466  }
467  Na[0] += sC_[i]*c[i];
468  Na[1] += sH_[i]*c[i];
469  Na[2] += sO_[i]*c[i];
470  if (sC_[i]>nbCLarge_ || this->chemistry_.Y()[i].member() == "O2")
471  {
472  Nal[0] += sC_[i]*c[i];
473  Nal[1] += sH_[i]*c[i];
474  Nal[2] += sO_[i]*c[i];
475  }
476  }
477 
478  // 2C(-CO2) + H(-H2O)/2 - z'C(-CO2)
479  // Progress equivalence ratio = ----------------------------------
480  // O(-CO2-H2O) - z' C(-CO2)
481  // where minus means that this species is not considered for the number
482  // of atoms and z' is the ratio of the number of O and C in the fuel(s)
483  phiProgress = (2*Na[0]+Na[1]/2-zprime_*Na[0])/(Na[2]-zprime_*Na[0]);
484 
485  // 2Cl + Hl/2
486  // Equivalence ratio for fuel decomposition = ----------
487  // Ol(+O2)
488  phiLarge = (2*Nal[0]+Nal[1]/2)/Nal[2];
489  }
490 
491  // Using the rAB matrix (numerator and denominator separated)
492  // compute the R value according to the search initiating set
493  scalarField Rvalue(this->nSpecie_,0.0);
494  label speciesNumber = 0;
495 
496  // Set all species to inactive and activate them according
497  // to rAB and initial set
498  for (label i=0; i<this->nSpecie_; i++)
499  {
500  this->activeSpecies_[i] = false;
501  }
502 
503  // Initialise the FIFOStack for search set
505 
506  const labelList& SIS(searchInitSet_);
507 
508  // If automaticSIS is on, the search initiating set is selected according to
509  // phiProgress and phiLarge
510  if (automaticSIS_)
511  {
512  if (phiLarge >= phiTol_ && phiProgress >= phiTol_)
513  {
514  // When phiLarge and phiProgress >= phiTol then
515  // CO, HO2 and fuel are in the SIS
516  Q.push(COId_);
517  speciesNumber++;
518  this->activeSpecies_[COId_] = true;
519  Rvalue[COId_] = 1.0;
520  Q.push(HO2Id_);
521  speciesNumber++;
522  this->activeSpecies_[HO2Id_] = true;
523  Rvalue[HO2Id_] = 1.0;
524  forAll(fuelSpeciesID_,i)
525  {
526  Q.push(fuelSpeciesID_[i]);
527  speciesNumber++;
528  this->activeSpecies_[fuelSpeciesID_[i]] = true;
529  Rvalue[fuelSpeciesID_[i]] = 1.0;
530  }
531 
532  }
533  else if (phiLarge < phiTol_ && phiProgress >= phiTol_)
534  {
535  // When phiLarge < phiTol and phiProgress >= phiTol then
536  // CO, HO2 are in the SIS
537  Q.push(COId_);
538  speciesNumber++;
539  this->activeSpecies_[COId_] = true;
540  Rvalue[COId_] = 1.0;
541  Q.push(HO2Id_);
542  speciesNumber++;
543  this->activeSpecies_[HO2Id_] = true;
544  Rvalue[HO2Id_] = 1.0;
545 
546  if (forceFuelInclusion_)
547  {
548  forAll(fuelSpeciesID_,i)
549  {
550  Q.push(fuelSpeciesID_[i]);
551  speciesNumber++;
552  this->activeSpecies_[fuelSpeciesID_[i]] = true;
553  Rvalue[fuelSpeciesID_[i]] = 1.0;
554  }
555  }
556  }
557  else
558  {
559  // When phiLarge and phiProgress< phiTol then
560  // CO2, H2O are in the SIS
561  Q.push(CO2Id_);
562  speciesNumber++;
563  this->activeSpecies_[CO2Id_] = true;
564  Rvalue[CO2Id_] = 1.0;
565 
566  Q.push(H2OId_);
567  speciesNumber++;
568  this->activeSpecies_[H2OId_] = true;
569  Rvalue[H2OId_] = 1.0;
570  if (forceFuelInclusion_)
571  {
572  forAll(fuelSpeciesID_,i)
573  {
574  Q.push(fuelSpeciesID_[i]);
575  speciesNumber++;
576  this->activeSpecies_[fuelSpeciesID_[i]] = true;
577  Rvalue[fuelSpeciesID_[i]] = 1.0;
578  }
579  }
580  }
581 
582  if (T>NOxThreshold_ && NOId_!=-1)
583  {
584  Q.push(NOId_);
585  speciesNumber++;
586  this->activeSpecies_[NOId_] = true;
587  Rvalue[NOId_] = 1.0;
588  }
589  }
590  else // No automaticSIS => all species of the SIS are added
591  {
592  for (label i=0; i<SIS.size(); i++)
593  {
594  label q = SIS[i];
595  this->activeSpecies_[q] = true;
596  speciesNumber++;
597  Q.push(q);
598  Rvalue[q] = 1.0;
599  }
600  }
601 
602  // Execute the main loop for R-value
603  while (!Q.empty())
604  {
605  label u = Q.pop();
606  scalar Den = max(PA[u],CA[u]);
607  if (Den != 0)
608  {
609  for (label v=0; v<NbrABInit[u]; v++)
610  {
611  label otherSpec = rABOtherSpec(u, v);
612  scalar rAB = mag(rABNum(u, v))/Den;
613  if (rAB > 1)
614  {
615  rAB = 1;
616  }
617  // The direct link is weaker than the user-defined tolerance
618  if (rAB >= this->tolerance())
619  {
620  scalar Rtemp = Rvalue[u]*rAB;
621  // a link analysed previously is stronger
622  // the (composed) link is stronger than the user-defined
623  // tolerance
624  if ((Rvalue[otherSpec]<Rtemp) && (Rtemp>=this->tolerance()))
625  {
626  Q.push(otherSpec);
627  Rvalue[otherSpec] = Rtemp;
628  if (!this->activeSpecies_[otherSpec])
629  {
630  this->activeSpecies_[otherSpec] = true;
631  speciesNumber++;
632  }
633  }
634  }
635  }
636  }
637  }
638 
639  // Put a flag on the reactions containing at least one removed species
640  forAll(this->chemistry_.reactions(), i)
641  {
642  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
643  this->chemistry_.reactionsDisabled()[i] = false;
644 
645  forAll(R.lhs(), s)
646  {
647  label ss = R.lhs()[s].index;
648  if (!this->activeSpecies_[ss])
649  {
650  // Flag the reaction to disable it
651  this->chemistry_.reactionsDisabled()[i] = true;
652  break;
653  }
654  }
655  if (!this->chemistry_.reactionsDisabled()[i])
656  {
657  forAll(R.rhs(), s)
658  {
659  label ss = R.rhs()[s].index;
660  if (!this->activeSpecies_[ss])
661  {
662  // Flag the reaction to disable it
663  this->chemistry_.reactionsDisabled()[i] = true;
664  break;
665  }
666  }
667  }
668  }
669 
670  this->NsSimp_ = speciesNumber;
671  scalarField& simplifiedC(this->chemistry_.simplifiedC());
672  simplifiedC.setSize(this->NsSimp_+2);
673  DynamicList<label>& s2c(this->chemistry_.simplifiedToCompleteIndex());
674  s2c.setSize(this->NsSimp_);
675  Field<label>& c2s(this->chemistry_.completeToSimplifiedIndex());
676 
677  label j = 0;
678  for (label i=0; i<this->nSpecie_; i++)
679  {
680  if (this->activeSpecies_[i])
681  {
682  s2c[j] = i;
683  simplifiedC[j] = c[i];
684  c2s[i] = j++;
685  if (!this->chemistry_.active(i))
686  {
687  this->chemistry_.setActive(i);
688  }
689  }
690  else
691  {
692  c2s[i] = -1;
693  }
694  }
695 
696  simplifiedC[this->NsSimp_] = T;
697  simplifiedC[this->NsSimp_+1] = p;
698  this->chemistry_.setNsDAC(this->NsSimp_);
699 
700  // Change temporary Ns in chemistryModel
701  // to make the function nEqns working
702  this->chemistry_.setNSpecie(this->NsSimp_);
703 }
704 
705 
706 // ************************************************************************* //
A FIFO stack based on a singly-linked list.
Definition: FIFOStack.H:51
Extends standardChemistryModel by adding the TDAC method.
A HashTable with keys but without contents.
Definition: HashSet.H:59
dictionary dict
#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
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:124
error FatalError
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:323
virtual void reduceMechanism(const scalar p, const scalar T, const scalarField &c, const label li)
Reduce the mechanism.
Definition: DAC.C:236
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: HashTable.H:59
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
const List< specieCoeffs > & lhs() const
Return the components of the left hand side.
Definition: reactionI.H:36
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:55
void omega(const scalar p, const scalar T, const scalarField &c, const label li, scalarField &dcdt) const
Net reaction rate for individual species.
Definition: Reaction.C:290
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:53
const dimensionedScalar c1
First radiation constant: default SI units: [W/m^2].
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))
void setSize(const label)
Alter the addressed list size.
Definition: DynamicListI.H:175
A class for handling words, derived from string.
Definition: word.H:59
const word & name() const
Return the name of the element.
label nAtoms() const
Return the number of atoms of this element in the specie.
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
basicChemistryModel & chemistry
const volScalarField & T
virtual ~DAC()
Destructor.
Definition: DAC.C:228
void setSize(const label)
Reset size of List.
Definition: List.C:281
DAC(const IOdictionary &dict, TDACChemistryModel< ThermoType > &chemistry)
Construct from components.
Definition: DAC.C:32
const List< specieCoeffs > & rhs() const
Return the components of the right hand side.
Definition: reactionI.H:42
#define R(A, B, C, D, E, F, K, M)
void push(const T &a)
Push an element onto the stack.
Definition: FIFOStack.H:84
messageStream Info
const multiComponentMixture< ThermoType > & mixture() const
Return reference to the mixture.
dimensioned< scalar > mag(const dimensioned< Type > &)
T pop()
Pop the bottom element off the stack.
Definition: FIFOStack.H:90
An abstract class for methods of chemical mechanism reduction.
volScalarField & p
bool found