Cantera  2.1.2
ReactorNet.h
Go to the documentation of this file.
1 /**
2  * @file ReactorNet.h
3  */
4 // Copyright 2004 California Institute of Technology
5 
6 #ifndef CT_REACTORNET_H
7 #define CT_REACTORNET_H
8 
9 #include "Reactor.h"
12 #include "cantera/base/Array.h"
13 
14 namespace Cantera
15 {
16 
17 //! A class representing a network of connected reactors.
18 /*!
19  * This class is used to integrate the time-dependent governing equations for
20  * a network of reactors (Reactor, ConstPressureReactor) connected by various
21  * means, e.g. Wall, MassFlowController, Valve, PressureController.
22  */
24 {
25 public:
26  ReactorNet();
27  virtual ~ReactorNet();
28 
29  /** @name Methods to set up a simulation. */
30  //@{
31 
32  /**
33  * Set initial time. Default = 0.0 s. Restarts integration
34  * from this time using the current mixture state as the
35  * initial condition.
36  */
37  void setInitialTime(doublereal time) {
38  m_time = time;
39  m_init = false;
40  }
41 
42  //! Set the maximum time step.
43  void setMaxTimeStep(double maxstep) {
44  m_maxstep = maxstep;
45  m_init = false;
46  }
47 
48  //! Set the maximum number of error test failures permitted by the CVODES
49  //! integrator in a single time step.
50  void setMaxErrTestFails(int nmax) {
51  m_maxErrTestFails = nmax;
52  m_init = false;
53  }
54 
55  //! Set the relative and absolute tolerances for the integrator.
56  void setTolerances(doublereal rtol, doublereal atol) {
57  if (rtol >= 0.0) {
58  m_rtol = rtol;
59  }
60  if (atol >= 0.0) {
61  m_atols = atol;
62  }
63  m_init = false;
64  }
65 
66  //! Set the relative and absolute tolerances for integrating the
67  //! sensitivity equations.
68  void setSensitivityTolerances(doublereal rtol, doublereal atol) {
69  if (rtol >= 0.0) {
70  m_rtolsens = rtol;
71  }
72  if (atol >= 0.0) {
73  m_atolsens = atol;
74  }
75  m_init = false;
76  }
77 
78  //! Current value of the simulation time.
79  doublereal time() {
80  return m_time;
81  }
82 
83  //! Relative tolerance.
84  doublereal rtol() {
85  return m_rtol;
86  }
87 
88  //! Absolute integration tolerance
89  doublereal atol() {
90  return m_atols;
91  }
92 
93  //! Relative sensitivity tolerance
94  doublereal rtolSensitivity() const {
95  return m_rtolsens;
96  }
97 
98  //! Absolute sensitivity tolerance
99  doublereal atolSensitivity() const {
100  return m_atolsens;
101  }
102 
103  /**
104  * Advance the state of all reactors in time. Take as many internal
105  * timesteps as necessary to reach *time*.
106  * @param time Time to advance to (s).
107  */
108  void advance(doublereal time);
109 
110  //! Advance the state of all reactors in time. Take a single timestep
111  //! toward *time*.
112  double step(doublereal time);
113 
114  //@}
115 
116  //! Add the reactor *r* to this reactor network.
117  void addReactor(ReactorBase* r, bool iown = false);
118 
119  //! Return a reference to the *n*-th reactor in this network. The reactor
120  //! indices are determined by the order in which the reactors were added
121  //! to the reactor network.
123  return *m_r[n];
124  }
125 
126  //! Returns `true` if verbose logging output is enabled.
127  bool verbose() const {
128  return m_verbose;
129  }
130 
131  //! Enable or disable verbose logging while setting up and integrating the
132  //! reactor network.
133  void setVerbose(bool v = true) {
134  m_verbose = v;
135  }
136 
137  //! Return a reference to the integrator.
139  return *m_integ;
140  }
141 
142  //! Update the state of all the reactors in the network to correspond to
143  //! the values in the solution vector *y*.
144  void updateState(doublereal* y);
145 
146  //! Return the sensitivity of the *k*-th solution component with
147  //! respect to the *p*-th sensitivity parameter.
148  double sensitivity(size_t k, size_t p) {
149  if (!m_init) {
150  initialize();
151  }
152  return m_integ->sensitivity(k, m_sensIndex[p])/m_integ->solution(k);
153  }
154 
155  //! Return the sensitivity of the species named *species* with respect to
156  //! the *p*-th sensitivity parameter.
157  double sensitivity(const std::string& species, size_t p, int reactor=0) {
158  size_t k = globalComponentIndex(species, reactor);
159  return sensitivity(k, p);
160  }
161 
162  //! Evaluate the Jacobian matrix for the reactor network.
163  /*!
164  * @param[in] t Time at which to evaluate the Jacobian
165  * @param[in] y Global state vector at time *t*
166  * @param[out] ydot Time derivative of the state vector evaluated at *t*.
167  * @param[in] p sensitivity parameter vector (unused?)
168  * @param[out] j Jacobian matrix, size neq() by neq().
169  */
170  void evalJacobian(doublereal t, doublereal* y,
171  doublereal* ydot, doublereal* p, Array2D* j);
172 
173  // overloaded methods of class FuncEval
174  virtual size_t neq() {
175  return m_nv;
176  }
177  virtual void eval(doublereal t, doublereal* y,
178  doublereal* ydot, doublereal* p);
179  virtual void getInitialConditions(doublereal t0, size_t leny,
180  doublereal* y);
181  virtual size_t nparams() {
182  return m_ntotpar;
183  }
184 
185  //! Return the index corresponding to the species named *species* in the
186  //! reactor with index *reactor* in the global state vector for the
187  //! reactor network.
188  size_t globalComponentIndex(const std::string& species, size_t reactor=0);
189 
190  //! Used by Reactor and Wall objects to register the addition of
191  //! sensitivity reactions so that the ReactorNet can keep track of the
192  //! order in which sensitivity parameters are added.
193  void registerSensitivityReaction(void* reactor, size_t reactionIndex,
194  const std::string& name, int leftright=0);
195 
196  //! The name of the p-th sensitivity parameter added to this ReactorNet.
197  const std::string& sensitivityParameterName(size_t p) {
198  return m_paramNames.at(p);
199  }
200 
201  void connect(size_t i, size_t j) {
202  m_connect[j*m_nr + i] = 1;
203  m_connect[i*m_nr + j] = 1;
204  }
205 
206  bool connected(size_t i, size_t j) {
207  return (m_connect[m_nr*i + j] == 1);
208  }
209 
210 protected:
211  /**
212  * Initialize the reactor network. Called automatically the first time
213  * advance or step is called.
214  */
215  void initialize();
216 
217  std::vector<ReactorBase*> m_r;
218  std::vector<Reactor*> m_reactors;
219  size_t m_nr;
220  size_t m_nreactors;
221  Integrator* m_integ;
222  doublereal m_time;
223  bool m_init;
224  size_t m_nv;
225  std::vector<size_t> m_size;
226  vector_fp m_atol;
227  doublereal m_rtol, m_rtolsens;
228  doublereal m_atols, m_atolsens;
229  doublereal m_maxstep;
230  int m_maxErrTestFails;
231  bool m_verbose;
232  size_t m_ntotpar;
233  std::vector<size_t> m_nparams;
234 
235  //! Names corresponding to each sensitivity parameter
236  std::vector<std::string> m_paramNames;
237 
238  //! Structure used to determine the order of sensitivity parameters
239  //! m_sensOrder[Reactor or Wall, leftright][reaction number] = parameter index
240  std::map<std::pair<void*, int>, std::map<size_t, size_t> > m_sensOrder;
241 
242  //! Mapping from the order in which sensitivity parameters were added to
243  //! the ReactorNet to the order in which they occur in the integrator
244  //! output.
245  std::vector<size_t> m_sensIndex;
246 
247  vector_int m_connect;
248  vector_fp m_ydot;
249 
250  std::vector<bool> m_iown;
251 };
252 }
253 
254 #endif
std::vector< size_t > m_sensIndex
Mapping from the order in which sensitivity parameters were added to the ReactorNet to the order in w...
Definition: ReactorNet.h:245
ReactorBase & reactor(int n)
Return a reference to the n-th reactor in this network.
Definition: ReactorNet.h:122
double step(doublereal time)
Advance the state of all reactors in time.
Definition: ReactorNet.cpp:161
A class representing a network of connected reactors.
Definition: ReactorNet.h:23
std::vector< std::string > m_paramNames
Names corresponding to each sensitivity parameter.
Definition: ReactorNet.h:236
void updateState(doublereal *y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
Definition: ReactorNet.cpp:238
double sensitivity(const std::string &species, size_t p, int reactor=0)
Return the sensitivity of the species named species with respect to the p-th sensitivity parameter...
Definition: ReactorNet.h:157
void setTolerances(doublereal rtol, doublereal atol)
Set the relative and absolute tolerances for the integrator.
Definition: ReactorNet.h:56
doublereal atolSensitivity() const
Absolute sensitivity tolerance.
Definition: ReactorNet.h:99
doublereal rtolSensitivity() const
Relative sensitivity tolerance.
Definition: ReactorNet.h:94
virtual size_t neq()
Number of equations.
Definition: ReactorNet.h:174
size_t globalComponentIndex(const std::string &species, size_t reactor=0)
Return the index corresponding to the species named species in the reactor with index reactor in the ...
Definition: ReactorNet.cpp:257
void initialize()
Initialize the reactor network.
Definition: ReactorNet.cpp:47
void evalJacobian(doublereal t, doublereal *y, doublereal *ydot, doublereal *p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
Definition: ReactorNet.cpp:204
Integrator & integrator()
Return a reference to the integrator.
Definition: ReactorNet.h:138
void setMaxTimeStep(double maxstep)
Set the maximum time step.
Definition: ReactorNet.h:43
A class for 2D arrays stored in column-major (Fortran-compatible) form.
Definition: Array.h:33
virtual void getInitialConditions(doublereal t0, size_t leny, doublereal *y)
Fill the solution vector with the initial conditions at initial time t0.
Definition: ReactorNet.cpp:247
Header file for class Cantera::Array2D.
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:167
void setVerbose(bool v=true)
Enable or disable verbose logging while setting up and integrating the reactor network.
Definition: ReactorNet.h:133
doublereal atol()
Absolute integration tolerance.
Definition: ReactorNet.h:89
void addReactor(ReactorBase *r, bool iown=false)
Add the reactor r to this reactor network.
Definition: ReactorNet.cpp:174
void setMaxErrTestFails(int nmax)
Set the maximum number of error test failures permitted by the CVODES integrator in a single time ste...
Definition: ReactorNet.h:50
Abstract base class for ODE system integrators.
Definition: Integrator.h:53
Base class for stirred reactors.
Definition: ReactorBase.h:30
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the right-hand-side function.
Definition: ReactorNet.cpp:188
void setSensitivityTolerances(doublereal rtol, doublereal atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
Definition: ReactorNet.h:68
virtual size_t nparams()
Number of sensitivity parameters.
Definition: ReactorNet.h:181
void registerSensitivityReaction(void *reactor, size_t reactionIndex, const std::string &name, int leftright=0)
Used by Reactor and Wall objects to register the addition of sensitivity reactions so that the Reacto...
Definition: ReactorNet.cpp:270
doublereal rtol()
Relative tolerance.
Definition: ReactorNet.h:84
bool verbose() const
Returns true if verbose logging output is enabled.
Definition: ReactorNet.h:127
std::vector< double > vector_fp
Turn on the use of stl vectors for the basic array type within cantera Vector of doubles.
Definition: ct_defs.h:165
virtual doublereal & solution(size_t k)
The current value of the solution of equation k.
Definition: Integrator.h:137
doublereal time()
Current value of the simulation time.
Definition: ReactorNet.h:79
std::map< std::pair< void *, int >, std::map< size_t, size_t > > m_sensOrder
Structure used to determine the order of sensitivity parameters m_sensOrder[Reactor or Wall...
Definition: ReactorNet.h:240
void setInitialTime(doublereal time)
Set initial time.
Definition: ReactorNet.h:37
Virtual base class for ODE right-hand-side function evaluators.
Definition: FuncEval.h:23
double sensitivity(size_t k, size_t p)
Return the sensitivity of the k-th solution component with respect to the p-th sensitivity parameter...
Definition: ReactorNet.h:148
const std::string & sensitivityParameterName(size_t p)
The name of the p-th sensitivity parameter added to this ReactorNet.
Definition: ReactorNet.h:197
void advance(doublereal time)
Advance the state of all reactors in time.
Definition: ReactorNet.cpp:148