PFA.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 "PFA.H"
27 
28 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
29 
30 template<class ThermoType>
32 (
33  const IOdictionary& dict,
35 )
36 :
38  searchInitSet_()
39 {
40  const wordHashSet initSet(this->coeffsDict_.lookup("initialSet"));
41  forAllConstIter(wordHashSet, initSet, iter)
42  {
43  searchInitSet_.append(chemistry.mixture().species()[iter.key()]);
44  }
45 }
46 
47 
48 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
49 
50 template<class ThermoType>
52 {}
53 
54 
55 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
56 
57 template<class ThermoType>
59 (
60  const scalar p,
61  const scalar T,
62  const scalarField& c,
63  const label li
64 )
65 {
66  scalarField& completeC(this->chemistry_.completeC());
67  scalarField c1(this->chemistry_.nEqns(), 0.0);
68 
69  for (label i=0; i<this->nSpecie_; i++)
70  {
71  c1[i] = c[i];
72  completeC[i] = c[i];
73  }
74 
75  c1[this->nSpecie_] = T;
76  c1[this->nSpecie_+1] = p;
77 
78  // Compute the rAB matrix
79  RectangularMatrix<scalar> PAB(this->nSpecie_,this->nSpecie_,0.0);
80  RectangularMatrix<scalar> CAB(this->nSpecie_,this->nSpecie_,0.0);
81  scalarField PA(this->nSpecie_,0.0);
82  scalarField CA(this->nSpecie_,0.0);
83 
84  // Number of initialised rAB for each lines
85  Field<label> NbrABInit(this->nSpecie_,0);
86  // Position of the initialised rAB, -1 when not initialised
87  RectangularMatrix<label> rABPos(this->nSpecie_, this->nSpecie_, -1);
88  // Index of the other species involved in the rABNum
89  RectangularMatrix<label> rABOtherSpec(this->nSpecie_, this->nSpecie_, -1);
90 
91  forAll(this->chemistry_.reactions(), i)
92  {
93  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
94 
95  // for each reaction compute omegai
96  scalar omegaf, omegar;
97  const scalar omegai = R.omega(p, T, c1, li, omegaf, omegar);
98 
99  // then for each pair of species composing this reaction,
100  // compute the rAB matrix (separate the numerator and
101  // denominator)
102 
103  DynamicList<scalar> wA(R.lhs().size()+R.rhs().size());
104  DynamicList<label> wAID(R.lhs().size()+R.rhs().size());
105 
106  forAll(R.lhs(), s)// compute rAB for all species in the left hand side
107  {
108  label ss = R.lhs()[s].index;
109  scalar sl = -R.lhs()[s].stoichCoeff; // vAi = v''-v' => here -v'
110  bool found(false);
111  forAll(wAID, id)
112  {
113  if (ss==wAID[id])
114  {
115  wA[id] += sl*omegai;
116  found = true;
117  break;
118  }
119  }
120  if (!found)
121  {
122  wA.append(sl*omegai);
123  wAID.append(ss);
124  }
125  }
126  forAll(R.rhs(), s) // compute rAB for all species in the right hand side
127  {
128  label ss = R.rhs()[s].index;
129  scalar sl = R.rhs()[s].stoichCoeff; // vAi = v''-v' => here v''
130  bool found(false);
131  forAll(wAID, id)
132  {
133  if (ss==wAID[id])
134  {
135  wA[id] += sl*omegai;
136  found = true;
137  break;
138  }
139  }
140  if (!found)
141  {
142  wA.append(sl*omegai);
143  wAID.append(ss);
144  }
145  }
146 
147  wAID.shrink();
148  forAll(wAID, id)
149  {
150  label curID = wAID[id];
151  scalar curwA = wA[id];
152  List<bool> deltaBi(this->nSpecie_, false);
153  FIFOStack<label> usedIndex;
154  forAll(R.lhs(),j)
155  {
156  label sj = R.lhs()[j].index;
157  usedIndex.push(sj);
158  deltaBi[sj] = true;
159  }
160  forAll(R.rhs(),j)
161  {
162  label sj = R.rhs()[j].index;
163  usedIndex.push(sj);
164  deltaBi[sj] = true;
165  }
166 
167  deltaBi[curID] = false;
168 
169  while(!usedIndex.empty())
170  {
171  label curIndex = usedIndex.pop();
172 
173  if (deltaBi[curIndex])
174  {
175  deltaBi[curIndex] = false;
176  if (rABPos(curID, curIndex)==-1)
177  {
178  rABPos(curID, curIndex) = NbrABInit[curID];
179  rABOtherSpec(curID, NbrABInit[curID]) = curIndex;
180  if (curwA > 0.0)
181  {
182  PAB(curID, NbrABInit[curID]) = curwA;
183  }
184  else
185  {
186  CAB(curID, NbrABInit[curID]) = -curwA;
187  }
188  NbrABInit[curID]++;
189  }
190  else
191  {
192  if (curwA > 0.0)
193  {
194  PAB(curID, rABPos(curID, curIndex)) += curwA;
195  }
196  else
197  {
198  CAB(curID, rABPos(curID, curIndex)) += -curwA;
199  }
200  }
201  }
202  }
203  // Now that every species of the reactions has been visited, we can
204  // compute the production and consumption rate. It avoids getting
205  // wrong results when species are present in both lhs and rhs
206  if (curwA > 0.0)
207  {
208  if (PA[curID] == 0.0)
209  {
210  PA[curID] = curwA;
211  }
212  else
213  {
214  PA[curID] += curwA;
215  }
216  }
217  else
218  {
219  if (CA[curID] == 0.0)
220  {
221  CA[curID] = -curwA;
222  }
223  else
224  {
225  CA[curID] += -curwA;
226  }
227  }
228  }
229  }
230  // rii = 0.0 by definition
231 
232  // Compute second generation link strength
233  // For all species A look at all rAri of the connected species ri and
234  // compute rriB with all the connected species of ri, B different of A. If
235  // a new species is connected, add it to the list of connected species. It
236  // is a connection of second generation and it will be aggregated in the
237  // final step to evaluate the total connection strength (or path flux).
238  // Compute rsecond=rAri*rriB with A!=ri!=B
239  RectangularMatrix<scalar> PAB2nd(this->nSpecie_,this->nSpecie_,0.0);
240  RectangularMatrix<scalar> CAB2nd(this->nSpecie_,this->nSpecie_,0.0);
241 
242  // Number of initialised rAB for each lines
243  Field<label> NbrABInit2nd(this->nSpecie_, 0);
244 
245  // Position of the initialised rAB, -1 when not initialised
246  RectangularMatrix<label> rABPos2nd(this->nSpecie_, this->nSpecie_, -1);
247 
248  // Index of the other species involved in the rABNum
249  RectangularMatrix<label> rABOtherSpec2nd
250  (
251  this->nSpecie_, this->nSpecie_, -1
252  );
253 
254  forAll(NbrABInit, A)
255  {
256  for (int i=0; i<NbrABInit[A]; i++)
257  {
258  label ri = rABOtherSpec(A, i);
259  scalar maxPACA = max(PA[ri],CA[ri]);
260  if (maxPACA > vSmall)
261  {
262  for (int j=0; j<NbrABInit[ri]; j++)
263  {
264  label B = rABOtherSpec(ri, j);
265  if (B != A) // if B!=A and also !=ri by definition
266  {
267  if (rABPos2nd(A, B)==-1)
268  {
269  rABPos2nd(A, B) = NbrABInit2nd[A]++;
270  rABOtherSpec2nd(A, rABPos2nd(A, B)) = B;
271  PAB2nd(A, rABPos2nd(A, B)) =
272  PAB(A, i)*PAB(ri, j)/maxPACA;
273  CAB2nd(A, rABPos2nd(A, B)) =
274  CAB(A, i)*CAB(ri, j)/maxPACA;
275  }
276  else
277  {
278  PAB2nd(A, rABPos2nd(A, B)) +=
279  PAB(A, i)*PAB(ri, j)/maxPACA;
280  CAB2nd(A, rABPos2nd(A, B)) +=
281  CAB(A, i)*CAB(ri, j)/maxPACA;
282  }
283  }
284  }
285  }
286  }
287  }
288 
289  // Using the rAB matrix (numerator and denominator separated)
290  label speciesNumber = 0;
291 
292  // set all species to inactive and activate them according
293  // to rAB and initial set
294  for (label i=0; i<this->nSpecie_; i++)
295  {
296  this->activeSpecies_[i] = false;
297  }
298 
299  // Initialise the FIFOStack for search set
300  const labelList& SIS(this->searchInitSet_);
302 
303  for (label i=0; i<SIS.size(); i++)
304  {
305  label q = SIS[i];
306  this->activeSpecies_[q] = true;
307  speciesNumber++;
308  Q.push(q);
309  }
310 
311  // Execute the main loop for R-value
312  while (!Q.empty())
313  {
314  label u = Q.pop();
315  scalar Den = max(PA[u],CA[u]);
316 
317  if (Den!=0.0)
318  {
319  // first generation
320  for (label v=0; v<NbrABInit[u]; v++)
321  {
322  label otherSpec = rABOtherSpec(u, v);
323  scalar rAB = (PAB(u, v)+CAB(u, v))/Den; // first generation
324  label id2nd = rABPos2nd(u, otherSpec);
325  if (id2nd !=-1)// if there is a second generation link
326  {
327  rAB += (PAB2nd(u, id2nd)+CAB2nd(u, id2nd))/Den;
328  }
329  // the link is stronger than the user-defined tolerance
330  if
331  (
332  rAB >= this->tolerance()
333  && !this->activeSpecies_[otherSpec]
334  )
335  {
336  Q.push(otherSpec);
337  this->activeSpecies_[otherSpec] = true;
338  speciesNumber++;
339  }
340 
341  }
342  // second generation link only (for those without first link)
343  for (label v=0; v<NbrABInit2nd[u]; v++)
344  {
345  label otherSpec = rABOtherSpec2nd(u, v);
346  scalar rAB = (PAB2nd(u, v)+CAB2nd(u, v))/Den;
347  // the link is stronger than the user-defined tolerance
348  if
349  (
350  rAB >= this->tolerance()
351  && !this->activeSpecies_[otherSpec]
352  )
353  {
354  Q.push(otherSpec);
355  this->activeSpecies_[otherSpec] = true;
356  speciesNumber++;
357  }
358  }
359  }
360  }
361 
362  // Put a flag on the reactions containing at least one removed species
363  forAll(this->chemistry_.reactions(), i)
364  {
365  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
366  this->chemistry_.reactionsDisabled()[i] = false;
367 
368  forAll(R.lhs(), s)
369  {
370  label ss = R.lhs()[s].index;
371  if (!this->activeSpecies_[ss])
372  {
373  this->chemistry_.reactionsDisabled()[i] = true;
374  break;
375  }
376  }
377  if (!this->chemistry_.reactionsDisabled()[i])
378  {
379  forAll(R.rhs(), s)
380  {
381  label ss = R.rhs()[s].index;
382  if (!this->activeSpecies_[ss])
383  {
384  this->chemistry_.reactionsDisabled()[i] = true;
385  break;
386  }
387  }
388  }
389  }
390 
391  this->NsSimp_ = speciesNumber;
392  scalarField& simplifiedC(this->chemistry_.simplifiedC());
393  simplifiedC.setSize(this->NsSimp_+2);
394  DynamicList<label>& s2c(this->chemistry_.simplifiedToCompleteIndex());
395  s2c.setSize(this->NsSimp_);
396  Field<label>& c2s(this->chemistry_.completeToSimplifiedIndex());
397 
398  label j = 0;
399  for (label i=0; i<this->nSpecie_; i++)
400  {
401  if (this->activeSpecies_[i])
402  {
403  s2c[j] = i;
404  simplifiedC[j] = c[i];
405  c2s[i] = j++;
406  if (!this->chemistry_.active(i))
407  {
408  this->chemistry_.setActive(i);
409  }
410  }
411  else
412  {
413  c2s[i] = -1;
414  }
415  }
416  simplifiedC[this->NsSimp_] = T;
417  simplifiedC[this->NsSimp_+1] = p;
418  this->chemistry_.setNsDAC(this->NsSimp_);
419  // change temporary Ns in chemistryModel
420  // to make the function nEqns working
421  this->chemistry_.setNSpecie(this->NsSimp_);
422 }
423 
424 
425 // ************************************************************************* //
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
dimensioned< Type > max(const dimensioned< Type > &, const dimensioned< Type > &)
void size(const label)
Override size to be inconsistent with allocated storage.
Definition: ListI.H:164
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
PFA(const IOdictionary &dict, TDACChemistryModel< ThermoType > &chemistry)
Construct from components.
Definition: PFA.C:32
forAllConstIter(PtrDictionary< phaseModel >, mixture.phases(), phase)
Definition: pEqn.H:29
basicChemistryModel & chemistry
virtual void reduceMechanism(const scalar p, const scalar T, const scalarField &c, const label li)
Reduce the mechanism.
Definition: PFA.C:59
const volScalarField & T
virtual ~PFA()
Destructor.
Definition: PFA.C:51
void setSize(const label)
Reset size of List.
Definition: List.C:281
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
const multiComponentMixture< ThermoType > & mixture() const
Return reference to the mixture.
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