DRG.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 -------------------------------------------------------------------------------
10 License
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 "DRG.H"
29 
30 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
31 
32 template<class CompType, class ThermoType>
34 (
35  const IOdictionary& dict,
37 )
38 :
40  searchInitSet_(this->coeffsDict_.subDict("initialSet").size())
41 {
42  label j=0;
43  dictionary initSet = this->coeffsDict_.subDict("initialSet");
44  for (label i=0; i<chemistry.nSpecie(); i++)
45  {
46  if (initSet.found(chemistry.Y()[i].member()))
47  {
48  searchInitSet_[j++] = i;
49  }
50  }
51  if (j<searchInitSet_.size())
52  {
54  << searchInitSet_.size()-j
55  << " species in the initial set is not in the mechanism "
56  << initSet
57  << exit(FatalError);
58  }
59 }
60 
61 
62 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
63 
64 template<class CompType, class ThermoType>
66 {}
67 
68 
69 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
70 
71 template<class CompType, class ThermoType>
73 (
74  const scalarField &c,
75  const scalar T,
76  const scalar p
77 )
78 {
79  scalarField c1(this->nSpecie_+2, Zero);
80 
81  for(label i=0; i<this->nSpecie_; i++)
82  {
83  c1[i] = c[i];
84  }
85 
86  c1[this->nSpecie_] = T;
87  c1[this->nSpecie_+1] = p;
88 
89  // Compute the rAB matrix
90  RectangularMatrix<scalar> rABNum(this->nSpecie_, this->nSpecie_, Zero);
91  scalarField rABDen(this->nSpecie_, Zero);
92 
93  // Number of initialized rAB for each lines
94  Field<label> NbrABInit(this->nSpecie_, Zero);
95 
96  // Position of the initialized rAB, -1 when not initialized
97  RectangularMatrix<label> rABPos(this->nSpecie_, this->nSpecie_, -1);
98 
99  // Index of the other species involved in the rABNum
100  RectangularMatrix<label> rABOtherSpec(this->nSpecie_, this->nSpecie_, -1);
101 
102  scalar pf, cf, pr, cr;
103  label lRef, rRef;
104  forAll(this->chemistry_.reactions(), i)
105  {
106  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
107  // For each reaction compute omegai
108  scalar omegai = this->chemistry_.omega
109  (
110  R, c1, T, p, pf, cf, lRef, pr, cr, rRef
111  );
112 
113 
114  // Then for each pair of species composing this reaction,
115  // compute the rAB matrix (separate the numerator and
116  // denominator)
117  DynamicList<scalar> wA(R.lhs().size()+R.rhs().size());
118  DynamicList<label> wAID(R.lhs().size()+R.rhs().size());
119 
120  forAll(R.lhs(), s)
121  {
122  label ss = R.lhs()[s].index;
123  scalar sl = -R.lhs()[s].stoichCoeff;
124  bool found(false);
125  forAll(wAID, id)
126  {
127  if (ss==wAID[id])
128  {
129  wA[id] += sl*omegai;
130  found = true;
131  break;
132  }
133  }
134  if (!found)
135  {
136  wA.append(sl*omegai);
137  wAID.append(ss);
138  }
139  }
140  forAll(R.rhs(), s)
141  {
142  label ss = R.rhs()[s].index;
143  scalar sl = R.rhs()[s].stoichCoeff;
144  bool found(false);
145  forAll(wAID, id)
146  {
147  if (ss==wAID[id])
148  {
149  wA[id] += sl*omegai;
150  found = true;
151  break;
152  }
153  }
154  if (!found)
155  {
156  wA.append(sl*omegai);
157  wAID.append(ss);
158  }
159  }
160 
161  // Now that all nuAi*wi are computed, without counting twice species
162  // present in both rhs and lhs, we can update the denominator and
163  // numerator for the rAB
164  wAID.shrink();
165  forAll(wAID, id)
166  {
167  label curID = wAID[id];
168 
169  // Absolute value of aggregated value
170  scalar curwA = ((wA[id]>=0) ? wA[id] : -wA[id]);
171 
172  List<bool> deltaBi(this->nSpecie_, false);
173  FIFOStack<label> usedIndex;
174  forAll(R.lhs(), j)
175  {
176  label sj = R.lhs()[j].index;
177  usedIndex.push(sj);
178  deltaBi[sj] = true;
179  }
180  forAll(R.rhs(), j)
181  {
182  label sj = R.rhs()[j].index;
183  usedIndex.push(sj);
184  deltaBi[sj] = true;
185  }
186 
187  // Disable for self reference (by definition rAA=0)
188  deltaBi[curID] = false;
189  while(!usedIndex.empty())
190  {
191  label curIndex = usedIndex.pop();
192 
193  if (deltaBi[curIndex])
194  {
195  // Disable to avoid counting it more than once
196  deltaBi[curIndex] = false;
197 
198  // Test if this rAB is not initialized
199  if (rABPos(curID, curIndex)==-1)
200  {
201  rABPos(curID, curIndex) = NbrABInit[curID];
202  NbrABInit[curID]++;
203  rABNum(curID, rABPos(curID, curIndex)) = curwA;
204  rABOtherSpec(curID, rABPos(curID, curIndex)) = curIndex;
205  }
206  else
207  {
208  rABNum(curID, rABPos(curID, curIndex)) += curwA;
209  }
210  }
211  }
212  if (rABDen[curID] == 0.0)
213  {
214  rABDen[curID] = curwA;
215  }
216  else
217  {
218  rABDen[curID] +=curwA;
219  }
220  }
221  }
222  // rii = 0.0 by definition
223 
224  label speciesNumber = 0;
225 
226  // Set all species to inactive and activate them according
227  // to rAB and initial set
228  for (label i=0; i<this->nSpecie_; i++)
229  {
230  this->activeSpecies_[i] = false;
231  }
232 
234 
235  // Initialize the list of active species with the search initiating set
236  // (SIS)
237  for (label i=0; i<searchInitSet_.size(); i++)
238  {
239  label q = searchInitSet_[i];
240  this->activeSpecies_[q] = true;
241  speciesNumber++;
242  Q.push(q);
243  }
244 
245  // Breadth first search with rAB
246  while (!Q.empty())
247  {
248  label u = Q.pop();
249  scalar Den = rABDen[u];
250 
251  if (Den > VSMALL)
252  {
253  for (label v=0; v<NbrABInit[u]; v++)
254  {
255  label otherSpec = rABOtherSpec(u, v);
256  scalar rAB = rABNum(u, v)/Den;
257 
258  if (rAB > 1)
259  {
260  rAB = 1;
261  }
262 
263  // Include B only if rAB is above the tolerance and if the
264  // species was not searched before
265  if
266  (
267  rAB >= this->tolerance()
268  && !this->activeSpecies_[otherSpec]
269  )
270  {
271  Q.push(otherSpec);
272  this->activeSpecies_[otherSpec] = true;
273  speciesNumber++;
274  }
275  }
276  }
277  }
278 
279  // Put a flag on the reactions containing at least one removed species
280  forAll(this->chemistry_.reactions(), i)
281  {
282  const Reaction<ThermoType>& R = this->chemistry_.reactions()[i];
283  this->chemistry_.reactionsDisabled()[i] = false;
284 
285  forAll(R.lhs(), s)
286  {
287  label ss = R.lhs()[s].index;
288 
289  // The species is inactive then the reaction is removed
290  if (!this->activeSpecies_[ss])
291  {
292  // Flag the reaction to disable it
293  this->chemistry_.reactionsDisabled()[i] = true;
294  break;
295  }
296  }
297 
298  // If the reaction has not been disabled yet
299  if (!this->chemistry_.reactionsDisabled()[i])
300  {
301  forAll(R.rhs(), s)
302  {
303  label ss = R.rhs()[s].index;
304  if (!this->activeSpecies_[ss])
305  {
306  this->chemistry_.reactionsDisabled()[i] = true;
307  break;
308  }
309  }
310  }
311  }
312 
313  this->NsSimp_ = speciesNumber;
314  this->chemistry_.simplifiedC().setSize(this->NsSimp_+2);
315  this->chemistry_.simplifiedToCompleteIndex().setSize(this->NsSimp_);
316 
317  label j = 0;
318  for (label i=0; i<this->nSpecie_; i++)
319  {
320  if (this->activeSpecies_[i])
321  {
322  this->chemistry_.simplifiedToCompleteIndex()[j] = i;
323  this->chemistry_.simplifiedC()[j] = c[i];
324  this->chemistry_.completeToSimplifiedIndex()[i] = j++;
325  if (!this->chemistry_.active(i))
326  {
327  this->chemistry_.setActive(i);
328  }
329  }
330  else
331  {
332  this->chemistry_.completeToSimplifiedIndex()[i] = -1;
333  }
334  }
335 
336  this->chemistry_.simplifiedC()[this->NsSimp_] = T;
337  this->chemistry_.simplifiedC()[this->NsSimp_+1] = p;
338  this->chemistry_.setNsDAC(this->NsSimp_);
339 
340  // Change temporary Ns in chemistryModel
341  // to make the function nEqns working
342  this->chemistry_.setNSpecie(this->NsSimp_);
343 }
344 
345 
346 // ************************************************************************* //
Foam::IOdictionary
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:54
DRG.H
Foam::chemistryReductionMethods::DRG::DRG
DRG(const IOdictionary &dict, TDACChemistryModel< CompType, ThermoType > &chemistry)
Construct from components.
Definition: DRG.C:34
p
volScalarField & p
Definition: createFieldRefs.H:8
Foam::chemistryReductionMethods::DRG::reduceMechanism
virtual void reduceMechanism(const scalarField &c, const scalar T, const scalar p)
Reduce the mechanism.
Definition: DRG.C:73
s
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))
Definition: gmvOutputSpray.H:25
Foam::Zero
static constexpr const zero Zero
Global zero (0)
Definition: zero.H:131
Foam::DynamicList< scalar >
Foam::dictionary::found
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
chemistry
BasicChemistryModel< psiReactionThermo > & chemistry
Definition: createFieldRefs.H:1
forAll
#define forAll(list, i)
Loop across all elements in list.
Definition: stdFoam.H:296
Foam::FIFOStack::push
void push(const T &element)
Push an element onto the back of the stack.
Definition: FIFOStack.H:78
R
#define R(A, B, C, D, E, F, K, M)
Foam::constant::physicoChemical::c1
const dimensionedScalar c1
First radiation constant: default SI units: [W/m2].
Foam::Field< scalar >
Foam::dictionary::subDict
const dictionary & subDict(const word &keyword, enum keyType::option matchOpt=keyType::REGEX) const
Find and return a sub-dictionary.
Definition: dictionary.C:460
Foam::RectangularMatrix< scalar >
dict
dictionary dict
Definition: searchingEngine.H:14
Foam::FatalError
error FatalError
Foam::dictionary
A list of keyword definitions, which are a keyword followed by a number of values (eg,...
Definition: dictionary.H:123
Foam::chemistryReductionMethods::DRG::~DRG
virtual ~DRG()
Definition: DRG.C:65
T
const volScalarField & T
Definition: createFieldRefs.H:2
Foam::exit
errorManipArg< error, int > exit(error &err, const int errNo=1)
Definition: errorManip.H:130
found
bool found
Definition: TABSMDCalcMethod2.H:32
Foam::FIFOStack::pop
T pop()
Pop the bottom element off the stack.
Definition: FIFOStack.H:90
FatalErrorInFunction
#define FatalErrorInFunction
Report an error message using Foam::FatalError.
Definition: error.H:453
Foam::List< bool >
Foam::constant::universal::c
const dimensionedScalar c
Speed of light in a vacuum.
Foam::TDACChemistryModel< CompType, ThermoType >
Foam::FIFOStack
A FIFO stack based on a singly-linked list.
Definition: FIFOStack.H:51
Foam::chemistryReductionMethod
An abstract class for methods of chemical mechanism reduction.
Definition: chemistryReductionMethod.H:56
Foam::Reaction
Simple extension of ReactionThermo to handle reaction kinetics in addition to the equilibrium thermod...
Definition: Reaction.H:59