DRGEP.C
Go to the documentation of this file.
1/*---------------------------------------------------------------------------*\
2 ========= |
3 \\ / F ield | OpenFOAM: The Open Source CFD Toolbox
4 \\ / O peration |
5 \\ / A nd | www.openfoam.com
6 \\/ M anipulation |
7-------------------------------------------------------------------------------
8 Copyright (C) 2016-2017 OpenFOAM Foundation
9-------------------------------------------------------------------------------
10License
11 This file is part of OpenFOAM.
12
13 OpenFOAM is free software: you can redistribute it and/or modify it
14 under the terms of the GNU General Public License as published by
15 the Free Software Foundation, either version 3 of the License, or
16 (at your option) any later version.
17
18 OpenFOAM is distributed in the hope that it will be useful, but WITHOUT
19 ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
20 FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
21 for more details.
22
23 You should have received a copy of the GNU General Public License
24 along with OpenFOAM. If not, see <http://www.gnu.org/licenses/>.
25
26\*---------------------------------------------------------------------------*/
27
28#include "DRGEP.H"
29#include "SortableListDRGEP.H"
30
31// * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
32
33template<class CompType, class ThermoType>
35(
36 const IOdictionary& dict,
38)
39:
40 chemistryReductionMethod<CompType, ThermoType>(dict, chemistry),
41 searchInitSet_(this->coeffsDict_.subDict("initialSet").size()),
42 sC_(this->nSpecie_, 0),
43 sH_(this->nSpecie_, 0),
44 sO_(this->nSpecie_, 0),
45 sN_(this->nSpecie_, 0),
46 NGroupBased_(50)
47{
48 label j = 0;
49 dictionary initSet = this->coeffsDict_.subDict("initialSet");
50 for (label i=0; i<chemistry.nSpecie(); ++i)
51 {
52 if (initSet.found(chemistry.Y()[i].member()))
53 {
54 searchInitSet_[j++] = i;
55 }
56 }
57 if (j<searchInitSet_.size())
58 {
60 << searchInitSet_.size() - j
61 << " species in the initial set is not in the mechanism "
62 << initSet
63 << exit(FatalError);
64 }
65
66 this->coeffsDict_.readIfPresent("NGroupBased", NGroupBased_);
67
68 const List<List<specieElement>>& specieComposition =
69 this->chemistry_.specieComp();
70
71 for (label i=0; i<this->nSpecie_; ++i)
72 {
73 const List<specieElement>& curSpecieComposition = specieComposition[i];
74
75 // for all elements in the current species
76 for (const specieElement& curElement : curSpecieComposition)
77 {
78 if (curElement.name() == "C")
79 {
80 sC_[i] = curElement.nAtoms();
81 }
82 else if (curElement.name() == "H")
83 {
84 sH_[i] = curElement.nAtoms();
85 }
86 else if (curElement.name() == "O")
87 {
88 sO_[i] = curElement.nAtoms();
89 }
90 else if (curElement.name() == "N")
91 {
92 sN_[i] = curElement.nAtoms();
93 }
94 else
95 {
96 Info<< "element not considered"<< endl;
97 }
98 }
99 }
100}
101
102
103// * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
104
105template<class CompType, class ThermoType>
107{}
108
109
110// * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
111
112template<class CompType, class ThermoType>
115(
116 const scalarField &c,
117 const scalar T,
118 const scalar p
119)
120{
121 scalarField& completeC(this->chemistry_.completeC());
122 scalarField c1(this->chemistry_.nEqns(), Zero);
123
124 for (label i=0; i<this->nSpecie_; ++i)
125 {
126 c1[i] = c[i];
127 completeC[i] = c[i];
128 }
129
130 c1[this->nSpecie_] = T;
131 c1[this->nSpecie_+1] = p;
132
133 // Compute the rAB matrix
134 RectangularMatrix<scalar> rABNum(this->nSpecie_, this->nSpecie_, Zero);
135 scalarField PA(this->nSpecie_, Zero);
136 scalarField CA(this->nSpecie_, Zero);
137
138 // Number of initialized rAB for each lines
139 Field<label> NbrABInit(this->nSpecie_, Zero);
140 // Position of the initialized rAB, -1 when not initialized
141 RectangularMatrix<label> rABPos(this->nSpecie_, this->nSpecie_, -1);
142 // Index of the other species involved in the rABNum
143 RectangularMatrix<label> rABOtherSpec(this->nSpecie_, this->nSpecie_, -1);
144
145 scalar pf, cf, pr, cr;
146 label lRef, rRef;
147 scalarField omegaV(this->chemistry_.reactions().size());
148 forAll(this->chemistry_.reactions(), i)
149 {
150 const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
151
152 // for each reaction compute omegai
153 scalar omegai = this->chemistry_.omega
154 (
155 R, c1, T, p, pf, cf, lRef, pr, cr, rRef
156 );
157 omegaV[i] = omegai;
158
159 // then for each pair of species composing this reaction,
160 // compute the rAB matrix (separate the numerator and
161 // denominator)
162
163 DynamicList<scalar> wA(R.lhs().size()+R.rhs().size());
164 DynamicList<label> wAID(R.lhs().size()+R.rhs().size());
165 forAll(R.lhs(), s)// compute rAB for all species in the left hand side
166 {
167 label ss = R.lhs()[s].index;
168 scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
169 List<bool> deltaBi(this->nSpecie_, false);
170 FIFOStack<label> usedIndex;
171 forAll(R.lhs(), j)
172 {
173 label sj = R.lhs()[j].index;
174 usedIndex.push(sj);
175 deltaBi[sj] = true;
176 }
177 forAll(R.rhs(), j)
178 {
179 label sj = R.rhs()[j].index;
180 usedIndex.push(sj);
181 deltaBi[sj] = true;
182 }
183
184 // Disable for self reference (by definition rAA=0)
185 deltaBi[ss] = false;
186
187 while (!usedIndex.empty())
188 {
189 label curIndex = usedIndex.pop();
190 if (deltaBi[curIndex])
191 {
192 // disable to avoid counting it more than once
193 deltaBi[curIndex] = false;
194 // test if this rAB is not initialized
195 if (rABPos(ss, curIndex) == -1)
196 {
197 rABPos(ss, curIndex) = NbrABInit[ss];
198 NbrABInit[ss]++;
199 rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
200 rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
201 }
202 else
203 {
204 rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
205 }
206 }
207 }
208 bool found(false);
209 forAll(wAID, id)
210 {
211 if (ss == wAID[id])
212 {
213 wA[id] += sl*omegai;
214 found = true;
215 }
216 }
217 if (!found)
218 {
219 wA.append(sl*omegai);
220 wAID.append(ss);
221 }
222 }
223
224 // Compute rAB for all species in the right hand side
225 forAll(R.rhs(), s)
226 {
227 label ss = R.rhs()[s].index;
228 scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
229 List<bool> deltaBi(this->nSpecie_, false);
230 FIFOStack<label> usedIndex;
231 forAll(R.lhs(), j)
232 {
233 label sj = R.lhs()[j].index;
234 usedIndex.push(sj);
235 deltaBi[sj] = true;
236 }
237 forAll(R.rhs(), j)
238 {
239 label sj = R.rhs()[j].index;
240 usedIndex.push(sj);
241 deltaBi[sj] = true;
242 }
243
244 // Disable for self reference (by definition rAA=0)
245 deltaBi[ss] = false;
246
247 while (!usedIndex.empty())
248 {
249 label curIndex = usedIndex.pop();
250 if (deltaBi[curIndex])
251 {
252 // disable to avoid counting it more than once
253 deltaBi[curIndex] = false;
254 // test if this rAB is not initialized
255 if (rABPos(ss, curIndex) == -1)
256 {
257 rABPos(ss, curIndex) = NbrABInit[ss];
258 NbrABInit[ss]++;
259 rABNum(ss, rABPos(ss, curIndex)) = sl*omegai;
260 rABOtherSpec(ss, rABPos(ss, curIndex)) = curIndex;
261 }
262 else
263 {
264 rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
265 }
266 }
267 }
268
269 bool found(false);
270 forAll(wAID, id)
271 {
272 if (ss == wAID[id])
273 {
274 wA[id] += sl*omegai;
275 found = true;
276 }
277 }
278 if (!found)
279 {
280 wA.append(sl*omegai);
281 wAID.append(ss);
282 }
283 }
284
285 wAID.shrink();
286 // Now that every species of the reactions has been visited, we can
287 // compute the production and consumption rate. This way, it avoids
288 // getting wrong results when species are present in both lhs and rhs
289 forAll(wAID, id)
290 {
291 if (wA[id] > 0.0)
292 {
293 if (PA[wAID[id]] == 0.0)
294 {
295 PA[wAID[id]] = wA[id];
296 }
297 else
298 {
299 PA[wAID[id]] += wA[id];
300 }
301 }
302 else
303 {
304 if (CA[wAID[id]] == 0.0)
305 {
306 CA[wAID[id]] = -wA[id];
307 }
308 else
309 {
310 CA[wAID[id]] += -wA[id];
311 }
312 }
313 }
314 }
315 // rii = 0.0 by definition
316
317 // Compute the production rate of each element Pa
318 label nElements = 4; // 4 main elements (C, H, O, N)
319 scalarList Pa(nElements, Zero);
320 scalarList Ca(nElements, Zero);
321
322 // for (label q=0; q<SIS.size(); ++q)
323 for (label i=0; i<this->nSpecie_; ++i)
324 {
325 Pa[0] += sC_[i]*max(0.0, (PA[i]-CA[i]));
326 Ca[0] += sC_[i]*max(0.0,-(PA[i]-CA[i]));
327 Pa[1] += sH_[i]*max(0.0, (PA[i]-CA[i]));
328 Ca[1] += sH_[i]*max(0.0,-(PA[i]-CA[i]));
329 Pa[2] += sO_[i]*max(0.0, (PA[i]-CA[i]));
330 Ca[2] += sO_[i]*max(0.0,-(PA[i]-CA[i]));
331 Pa[3] += sN_[i]*max(0.0, (PA[i]-CA[i]));
332 Ca[3] += sN_[i]*max(0.0,-(PA[i]-CA[i]));
333 }
334
335 // Using the rAB matrix (numerator and denominator separated)
336 // compute the R value according to the search initiating set
337 scalarField Rvalue(this->nSpecie_, Zero);
338 label speciesNumber = 0;
339 List<bool> disabledSpecies(this->nSpecie_, false);
340
341 // set all species to inactive and activate them according
342 // to rAB and initial set
343 for (label i=0; i<this->nSpecie_; ++i)
344 {
345 this->activeSpecies_[i] = false;
346 }
347 // Initialize the FIFOStack for search set
349 const labelList& SIS(this->searchInitSet_);
350 DynamicList<label> QStart(SIS.size());
351 DynamicList<scalar> alphaQ(SIS.size());
352
353 // Compute the alpha coefficient and initialize the R value of the species
354 // in the SIS
355 for (label i=0; i<SIS.size(); ++i)
356 {
357 label q = SIS[i];
358 // compute alpha
359 scalar alphaA(0.0);
360 // C
361 if (Pa[0] > VSMALL)
362 {
363 scalar alphaTmp = (sC_[q]*mag(PA[q]-CA[q])/Pa[0]);
364 if (alphaTmp > alphaA)
365 {
366 alphaA = alphaTmp;
367 }
368 }
369 // H
370 if (Pa[1] > VSMALL)
371 {
372 scalar alphaTmp = (sH_[q]*mag(PA[q]-CA[q])/Pa[1]);
373 if (alphaTmp > alphaA)
374 {
375 alphaA = alphaTmp;
376 }
377 }
378 // O
379 if (Pa[2] > VSMALL)
380 {
381 scalar alphaTmp = (sO_[q]*mag(PA[q]-CA[q])/Pa[2]);
382 if (alphaTmp > alphaA)
383 {
384 alphaA = alphaTmp;
385 }
386 }
387 // N
388 if (Pa[3] > VSMALL)
389 {
390 scalar alphaTmp = (sN_[q]*mag(PA[q]-CA[q])/Pa[3]);
391 if (alphaTmp > alphaA)
392 {
393 alphaA = alphaTmp;
394 }
395 }
396 if (alphaA > this->tolerance())
397 {
398 this->activeSpecies_[q] = true;
399 ++speciesNumber;
400 Q.push(q);
401 QStart.append(q);
402 alphaQ.append(1.0);
403 Rvalue[q] = 1.0;
404 }
405 else
406 {
407 Rvalue[q] = alphaA;
408 }
409 }
410
411 // if all species from the SIS has been removed
412 // force the use of the species with maximum Rvalue
413 if (Q.empty())
414 {
415 scalar Rmax=0.0;
416 label specID=-1;
417 for (const label sis : SIS)
418 {
419 if (Rvalue[sis] > Rmax)
420 {
421 Rmax = Rvalue[sis];
422 specID = sis;
423 }
424 }
425
426 Q.push(specID);
427 QStart.append(specID);
428 alphaQ.append(1.0);
429 ++speciesNumber;
430 Rvalue[specID] = 1.0;
431 this->activeSpecies_[specID] = true;
432 }
433
434 // Execute the main loop for R-value
435 while (!Q.empty())
436 {
437 label u = Q.pop();
438 scalar Den = max(PA[u], CA[u]);
439 if (Den > VSMALL)
440 {
441 for (label v=0; v<NbrABInit[u]; ++v)
442 {
443 label otherSpec = rABOtherSpec(u, v);
444 scalar rAB = mag(rABNum(u, v))/Den;
445 if (rAB > 1)
446 {
447 rAB = 1;
448 }
449
450 scalar Rtemp = Rvalue[u]*rAB;
451 // a link analysed previously is stronger
452 if (Rvalue[otherSpec] < Rtemp)
453 {
454 Rvalue[otherSpec] = Rtemp;
455 // the (composed) link is stronger than the tolerance
456 if (Rtemp >= this->tolerance())
457 {
458 Q.push(otherSpec);
459 if (!this->activeSpecies_[otherSpec])
460 {
461 this->activeSpecies_[otherSpec] = true;
462 ++speciesNumber;
463 }
464 }
465 }
466 }
467 }
468 }
469
470 // Group-based reduction
471 // number of species disabled in the first step
472 label NDisabledSpecies(this->nSpecie_ - speciesNumber);
473
474 // while the number of removed species is greater than NGroupBased, the rAB
475 // are reevaluated according to the group based definition for each loop the
476 // temporary disabled species (in the first reduction) are sorted to disable
477 // definitely the NGroupBased species with lower R then these R value a
478 // reevaluated taking into account these disabled species
479 while (NDisabledSpecies > NGroupBased_)
480 {
481 // keep track of disabled species using sortablelist to extract only
482 // NGroupBased lower R value
483 SortableListDRGEP<scalar> Rdisabled(NDisabledSpecies);
484 labelList Rindex(NDisabledSpecies);
485 label nD = 0;
486 forAll(disabledSpecies, i)
487 {
488 // if just disabled and not in a previous loop
489 if (!this->activeSpecies_[i] && !disabledSpecies[i])
490 {
491 // Note: non-reached species will be removed first (Rvalue=0)
492 Rdisabled[nD] = Rvalue[i];
493 Rindex[nD++] = i;
494 }
495 }
496 // sort the Rvalue to obtain the NGroupBased lower R value
497 Rdisabled.partialSort(NGroupBased_);
498 labelList tmpIndex(Rdisabled.indices());
499
500 // disable definitely NGroupBased species in this loop
501 for (label i=0; i<NGroupBased_; ++i)
502 {
503 disabledSpecies[Rindex[tmpIndex[i]]] = true;
504 }
505 NDisabledSpecies -= NGroupBased_;
506
507 // reevaluate the rAB according to the group-based definition rAB{S} [1]
508 // only update the numerator
509 forAll(NbrABInit, i)
510 {
511 for (label v=0; v<NbrABInit[i]; ++v)
512 {
513 rABNum(i, v) = 0.0;
514 }
515 }
516 forAll(this->chemistry_.reactions(), i)
517 {
518 const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
519
520 scalar omegai = omegaV[i];
521
522 forAll(R.lhs(), s)
523 {
524 label ss = R.lhs()[s].index;
525 scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
526 List<bool> deltaBi(this->nSpecie_, false);
527 bool alreadyDisabled(false);
528 FIFOStack<label> usedIndex;
529 forAll(R.lhs(), j)
530 {
531 label sj = R.lhs()[j].index;
532 usedIndex.push(sj);
533 deltaBi[sj] = true;
534 if (disabledSpecies[sj])
535 {
536 alreadyDisabled=true;
537 }
538 }
539 forAll(R.rhs(), j)
540 {
541 label sj = R.rhs()[j].index;
542 usedIndex.push(sj);
543 deltaBi[sj] = true;
544 if (disabledSpecies[sj])
545 {
546 alreadyDisabled=true;
547 }
548 }
549
550 deltaBi[ss] = false;
551
552 if (alreadyDisabled)
553 {
554 // if one of the species in this reaction is disabled, all
555 // species connected to species ss are modified
556 for (label v=0; v<NbrABInit[ss]; ++v)
557 {
558 rABNum(ss, v) += sl*omegai;
559 }
560 }
561 else
562 {
563 while(!usedIndex.empty())
564 {
565 label curIndex = usedIndex.pop();
566 if (deltaBi[curIndex])
567 {
568 // disable to avoid counting it more than once
569 deltaBi[curIndex] = false;
570 rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
571 }
572 }
573 }
574 }
575
576 forAll(R.rhs(), s)
577 {
578 label ss = R.rhs()[s].index;
579 scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
580 List<bool> deltaBi(this->nSpecie_, false);
581 bool alreadyDisabled(false);
582 FIFOStack<label> usedIndex;
583 forAll(R.lhs(), j)
584 {
585 label sj = R.lhs()[j].index;
586 usedIndex.push(sj);
587 deltaBi[sj] = true;
588 if (disabledSpecies[sj])
589 {
590 alreadyDisabled = true;
591 }
592 }
593 forAll(R.rhs(), j)
594 {
595 label sj = R.rhs()[j].index;
596 usedIndex.push(sj);
597 deltaBi[sj] = true;
598 if (disabledSpecies[sj])
599 {
600 alreadyDisabled = true;
601 }
602 }
603
604 deltaBi[ss] = false;
605
606 if (alreadyDisabled)
607 {
608 // if one of the species in this reaction is disabled, all
609 // species connected to species ss are modified
610 for (label v=0; v<NbrABInit[ss]; ++v)
611 {
612 rABNum(ss, v) += sl*omegai;
613 }
614 }
615 else
616 {
617 while(!usedIndex.empty())
618 {
619 label curIndex = usedIndex.pop();
620 if (deltaBi[curIndex])
621 {
622 deltaBi[curIndex] = false;
623 rABNum(ss, rABPos(ss, curIndex)) += sl*omegai;
624 }
625 }
626 }
627 }
628 }
629
630 for (const label q : QStart)
631 {
632 Q.push(q);
633 }
634
635 while (!Q.empty())
636 {
637 label u = Q.pop();
638 scalar Den = max(PA[u],CA[u]);
639 if (Den != 0.0)
640 {
641 for (label v=0; v<NbrABInit[u]; ++v)
642 {
643 label otherSpec = rABOtherSpec(u, v);
644 if (!disabledSpecies[otherSpec])
645 {
646 scalar rAB = mag(rABNum(u, v))/Den;
647 if (rAB > 1)
648 {
649 rAB = 1;
650 }
651
652 scalar Rtemp = Rvalue[u]*rAB;
653 // a link analysed previously is stronger
654 if (Rvalue[otherSpec] < Rtemp)
655 {
656 Rvalue[otherSpec] = Rtemp;
657 if (Rtemp >= this->tolerance())
658 {
659 Q.push(otherSpec);
660 if (!this->activeSpecies_[otherSpec])
661 {
662 this->activeSpecies_[otherSpec] = true;
663 ++speciesNumber;
664 NDisabledSpecies--;
665 }
666 }
667 }
668 }
669 }
670 }
671 }
672 }
673
674 // End of group-based reduction
675
676 // Put a flag on the reactions containing at least one removed species
677 forAll(this->chemistry_.reactions(), i)
678 {
679 const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
680 this->chemistry_.reactionsDisabled()[i] = false;
681 forAll(R.lhs(), s)
682 {
683 label ss = R.lhs()[s].index;
684 if (!this->activeSpecies_[ss])
685 {
686 this->chemistry_.reactionsDisabled()[i] = true;
687 break;
688 }
689 }
690 if (!this->chemistry_.reactionsDisabled()[i])
691 {
692 forAll(R.rhs(), s)
693 {
694 label ss = R.rhs()[s].index;
695 if (!this->activeSpecies_[ss])
696 {
697 this->chemistry_.reactionsDisabled()[i] = true;
698 break;
699 }
700 }
701 }
702 }
703
704 this->NsSimp_ = speciesNumber;
705 scalarField& simplifiedC(this->chemistry_.simplifiedC());
706 simplifiedC.setSize(this->NsSimp_+2);
707 DynamicList<label>& s2c(this->chemistry_.simplifiedToCompleteIndex());
708 s2c.setSize(this->NsSimp_);
709 Field<label>& c2s(this->chemistry_.completeToSimplifiedIndex());
710
711 label j = 0;
712 for (label i=0; i<this->nSpecie_; ++i)
713 {
714 if (this->activeSpecies_[i])
715 {
716 s2c[j] = i;
717 simplifiedC[j] = c[i];
718 c2s[i] = ++j;
719 if (!this->chemistry_.active(i))
720 {
721 this->chemistry_.setActive(i);
722 }
723 }
724 else
725 {
726 c2s[i] = -1;
727 }
728 }
729 simplifiedC[this->NsSimp_] = T;
730 simplifiedC[this->NsSimp_+1] = p;
731 this->chemistry_.setNsDAC(this->NsSimp_);
732 // change temporary Ns in chemistryModel
733 // to make the function nEqns working
734 this->chemistry_.setNSpecie(this->NsSimp_);
735}
736
737
738// ************************************************************************* //
#define R(A, B, C, D, E, F, K, M)
bool found
A 1D vector of objects of type <T> that resizes itself as necessary to accept the new objects.
Definition: DynamicList.H:72
DynamicList< T, SizeMin > & shrink()
Shrink the allocated space to the number of elements used.
Definition: DynamicListI.H:434
void setSize(const label n)
Same as resize()
Definition: DynamicList.H:221
void append(const T &val)
Copy append an element to the end of this list.
Definition: DynamicListI.H:503
A FIFO stack based on a singly-linked list.
Definition: FIFOStack.H:54
void push(const T &elem)
Push an element onto the back of the stack.
Definition: FIFOStack.H:78
T pop()
Pop the bottom element off the stack.
Definition: FIFOStack.H:90
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
A 1D array of objects of type <T>, where the size of the vector is known and used for subscript bound...
Definition: List.H:77
void setSize(const label n)
Alias for resize()
Definition: List.H:218
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:73
A templated (M x N) rectangular matrix of objects of <Type>, containing M*N elements,...
A list that is sorted upon construction or when explicitly requested with the sort() method.
const labelList & indices() const
Return the list of sorted indices. Updated every sort.
void partialSort(int M)
Partial sort the list (if changed after construction time)
Extends StandardChemistryModel by adding the TDAC method.
List< List< specieElement > > & specieComp()
void size(const label n)
Older name for setAddressableSize.
Definition: UList.H:114
An abstract class for methods of chemical mechanism reduction.
TDACChemistryModel< CompType, ThermoType > & chemistry_
const label nSpecie_
Number of species.
const dictionary coeffsDict_
Dictionary that store the algorithm data.
The DRGEP algorithm [1] is based on.
Definition: DRGEP.H:124
virtual ~DRGEP()
Destructor.
Definition: DRGEP.C:106
virtual void reduceMechanism(const scalarField &c, const scalar T, const scalar p)
Reduce the mechanism.
Definition: DRGEP.C:115
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:126
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
bool found(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Search for an entry (const access) with the given keyword.
Definition: dictionaryI.H:87
bool readIfPresent(const word &keyword, T &val, enum keyType::option matchOpt=keyType::REGEX) const
volScalarField & p
BasicChemistryModel< psiReactionThermo > & chemistry
const volScalarField & T
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
gmvFile<< "tracers "<< particles.size()<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().x()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().y()<< " ";}gmvFile<< nl;for(const passiveParticle &p :particles){ gmvFile<< p.position().z()<< " ";}gmvFile<< nl;forAll(lagrangianScalarNames, i){ word name=lagrangianScalarNames[i];IOField< scalar > s(IOobject(name, runTime.timeName(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
label max(const labelHashSet &set, label maxValue=labelMin)
Find the max value in labelHashSet, optionally limited by second argument.
Definition: hashSets.C:47
messageStream Info
Information stream (stdout output on master, null elsewhere)
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:372
dimensioned< typename typeOfMag< Type >::type > mag(const dimensioned< Type > &dt)
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
error FatalError
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
dictionary dict
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:333