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