chemistryReductionMethod.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) 2021-2022 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 
27 #include "chemistryModel.H"
28 
29 // * * * * * * * * * * * * * * * * Constructors * * * * * * * * * * * * * * //
30 
31 template<class ThermoType>
33 (
35 )
36 :
37  coeffsDict_(),
38  chemistry_(chemistry),
39  nSpecie_(chemistry.nSpecie()),
40  nActiveSpecies_(chemistry.nSpecie()),
41  reactionsDisabled_(chemistry.nReaction(), false),
42  activeSpecies_(chemistry.nSpecie(), true),
43  log_(false),
44  tolerance_(NaN),
45  sumnActiveSpecies_(0),
46  sumn_(0),
47  reduceMechCpuTime_(0)
48 {}
49 
50 
51 template<class ThermoType>
53 (
54  const Foam::IOdictionary& dict,
56 )
57 :
58  coeffsDict_(dict.subDict("reduction")),
59  chemistry_(chemistry),
60  nSpecie_(chemistry.nSpecie()),
61  nActiveSpecies_(chemistry.nSpecie()),
62  reactionsDisabled_(chemistry.nReaction(), false),
63  activeSpecies_(chemistry.nSpecie(), false),
64  log_(coeffsDict_.lookupOrDefault<Switch>("log", false)),
65  tolerance_(coeffsDict_.lookupOrDefault<scalar>("tolerance", 1e-4)),
66  sumnActiveSpecies_(0),
67  sumn_(0),
68  reduceMechCpuTime_(0)
69 {
70  if (log_)
71  {
72  cpuReduceFile_ = chemistry.logFile("cpu_reduce.out");
73  nActiveSpeciesFile_ = chemistry.logFile("nActiveSpecies.out");
74  }
75 }
76 
77 
78 // * * * * * * * * * * * * * * * * Destructor * * * * * * * * * * * * * * * //
79 
80 template<class ThermoType>
82 {}
83 
84 
85 // * * * * * * * * * * * * * * * Member Functions * * * * * * * * * * * * * //
86 
87 template<class ThermoType>
89 {
90  if (log_)
91  {
92  cpuTime_.cpuTimeIncrement();
93  }
94 }
95 
96 
97 template<class ThermoType>
99 (
100  List<label>& ctos,
101  DynamicList<label>& stoc
102 )
103 {
104  // Disable reactions containing removed species
105  forAll(chemistry_.reactions(), i)
106  {
107  const Reaction<ThermoType>& R = chemistry_.reactions()[i];
108  reactionsDisabled_[i] = false;
109 
110  forAll(R.lhs(), s)
111  {
112  label ss = R.lhs()[s].index;
113  if (!activeSpecies_[ss])
114  {
115  reactionsDisabled_[i] = true;
116  break;
117  }
118  }
119 
120  if (!reactionsDisabled_[i])
121  {
122  forAll(R.rhs(), s)
123  {
124  label ss = R.rhs()[s].index;
125  if (!activeSpecies_[ss])
126  {
127  reactionsDisabled_[i] = true;
128  break;
129  }
130  }
131  }
132  }
133 
134  // Set the total number of active species
135  nActiveSpecies_ = count(activeSpecies_, true);
136 
137  // Set the indexing arrays
138  stoc.setSize(nActiveSpecies_);
139  for (label i=0, j=0; i<nSpecie(); i++)
140  {
141  if (activeSpecies_[i])
142  {
143  stoc[j] = i;
144  ctos[i] = j++;
145  if (!chemistry_.active(i))
146  {
147  chemistry_.setActive(i);
148  }
149  }
150  else
151  {
152  ctos[i] = -1;
153  }
154  }
155 
156  // Change the number of species in the chemistry model
157  chemistry_.setNSpecie(nActiveSpecies_);
158 
159  if (log_)
160  {
161  sumnActiveSpecies_ += nActiveSpecies_;
162  sumn_++;
163  reduceMechCpuTime_ += cpuTime_.cpuTimeIncrement();
164  }
165 }
166 
167 
168 template<class ThermoType>
170 {
171  if (log_)
172  {
173  cpuReduceFile_()
174  << chemistry_.time().userTimeValue()
175  << " " << reduceMechCpuTime_ << endl;
176 
177  if (sumn_)
178  {
179  // Write average number of species
180  nActiveSpeciesFile_()
181  << chemistry_.time().userTimeValue()
182  << " " << sumnActiveSpecies_/sumn_ << endl;
183  }
184 
185  sumnActiveSpecies_ = 0;
186  sumn_ = 0;
187  reduceMechCpuTime_ = 0;
188  }
189 }
190 
191 
192 // ************************************************************************* //
#define forAll(list, i)
Loop across all elements in list.
Definition: UList.H:434
void setSize(const label)
Alter the addressed list size.
Definition: DynamicListI.H:175
IOdictionary is derived from dictionary and IOobject to give the dictionary automatic IO functionalit...
Definition: IOdictionary.H:57
Simple extension of ThermoType to handle reaction kinetics in addition to the equilibrium thermodynam...
Definition: Reaction.H:72
A simple wrapper around bool so that it can be read as a word: true/false, on/off,...
Definition: Switch.H:61
Extends base chemistry model by adding a thermo package, and ODE functions. Introduces chemistry equa...
virtual ~chemistryReductionMethod()
Destructor.
void initReduceMechanism()
Protected Member Functions.
void endReduceMechanism(List< label > &ctos, DynamicList< label > &stoc)
End reduction of the mechanism.
chemistryReductionMethod(chemistryModel< ThermoType > &chemistry)
Construct from components.
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.name(), cloud::prefix, mesh, IOobject::MUST_READ, IOobject::NO_WRITE))
static Type NaN()
Return a primitive with all components set to NaN.
const doubleScalar e
Definition: doubleScalar.H:105
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
Ostream & endl(Ostream &os)
Add newline and flush stream.
Definition: Ostream.H:251
static scalar R(const scalar a, const scalar x)
Definition: invIncGamma.C:102
label count(const ListType &l, typename ListType::const_reference x)
Count the number of occurrences of a value in a list.
const label nSpecie
dictionary dict
basicChemistryModel & chemistry