24ReactorNet::ReactorNet()
29ReactorNet::ReactorNet(shared_ptr<ReactorBase> reactor)
35ReactorNet::ReactorNet(vector<shared_ptr<ReactorBase>>& reactors)
38 for (
auto&
reactor : reactors) {
43ReactorNet::~ReactorNet()
91 throw CanteraError(
"ReactorNet::time",
"Time is not the independent variable"
92 " for this reactor network.");
100 throw CanteraError(
"ReactorNet::distance",
"Distance is not the independent"
101 " variable for this reactor network.");
108 debuglog(
"Initializing reactor network.\n", m_verbose);
109 if (m_reactors.empty()) {
111 "no reactors in network!");
114 for (
size_t n = 0; n < m_reactors.size(); n++) {
122 writelog(
"Reactor {:d}: {:d} variables.\n", n, nv);
125 if (r.
type() ==
"FlowReactor" && m_reactors.size() > 1) {
127 "FlowReactors must be used alone.");
131 m_ydot.resize(m_nv,0.0);
132 m_yest.resize(m_nv,0.0);
133 m_advancelimits.resize(m_nv,-1.0);
134 m_atol.resize(
neq());
135 fill(m_atol.begin(), m_atol.end(), m_atols);
136 m_integ->setTolerances(m_rtol,
neq(), m_atol.data());
137 m_integ->setSensitivityTolerances(m_rtolsens, m_atolsens);
138 if (!m_linearSolverType.empty()) {
139 m_integ->setLinearSolverType(m_linearSolverType);
142 m_integ->setPreconditioner(m_precon);
144 m_integ->initialize(
m_time, *
this);
149 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
159 debuglog(
"Re-initializing reactor network.\n", m_verbose);
160 m_integ->reinitialize(
m_time, *
this);
161 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
172 m_linearSolverType = linSolverType;
178 m_precon = preconditioner;
199 m_integ->integrate(
time);
232 double t =
time, delta;
233 double* y = m_integ->solution();
237 bool exceeded =
false;
239 for (
size_t j = 0; j < m_nv; j++) {
240 delta = abs(m_yest[j] - y[j]);
241 if ( (m_advancelimits[j] > 0.) && ( delta > m_advancelimits[j]) ) {
244 writelog(
" Limiting global state vector component {:d} (dt = {:9.4g}):"
245 "{:11.6g} > {:9.4g}\n",
246 j, t -
m_time, delta, m_advancelimits[j]);
278 vector<double> y(
neq());
282 solver.
solve(loglevel);
294 vector<double> y0(
neq());
295 vector<double> y1(
neq());
302 return std::dynamic_pointer_cast<EigenSparseJacobian>(solver.
linearSolver())->jacobian();
311 double* cvode_dky = m_integ->solution();
312 for (
size_t j = 0; j < m_nv; j++) {
313 yest[j] = cvode_dky[j];
319 for (
int n = 1; n <= k; n++) {
320 factor *= deltat / n;
321 cvode_dky = m_integ->derivative(
m_time, n);
322 for (
size_t j = 0; j < m_nv; j++) {
323 yest[j] += factor * cvode_dky[j];
331 return m_integ->lastOrder();
340 "To be removed after Cantera 3.2. Replaceable by reactor net "
341 "instantiation with contents.");
342 for (
auto current : m_reactors) {
343 if (current->isOde() != r.
isOde()) {
345 "Cannot mix Reactor types using both ODEs and DAEs ({} and {})",
346 current->type(), r.
type());
350 "Cannot mix Reactor types using time and space as independent variables"
351 "\n({} and {})", current->type(), r.
type());
356 m_reactors.push_back(&r);
362 m_integ->setLinearSolverType(
"DENSE");
369 auto r = std::dynamic_pointer_cast<Reactor>(
reactor);
372 "Reactor with type '{}' cannot be added to network.",
376 for (
auto current : m_reactors) {
377 if (current->isOde() != r->isOde()) {
379 "Cannot mix Reactor types using both ODEs and DAEs ({} and {})",
380 current->type(), r->type());
382 if (current->timeIsIndependent() != r->timeIsIndependent()) {
384 "Cannot mix Reactor types using time and space as independent variables"
385 "\n({} and {})", current->type(), r->type());
390 m_reactors.push_back(r.get());
396 m_integ->setLinearSolverType(
"DENSE");
406 for (
size_t i=0; i<r.
nWalls(); i++) {
409 if (w.left().type() ==
"Reservoir") {
412 if (w.right().type() ==
"Reservoir") {
417 for (
size_t i=0; i<r.
nInlets(); i++) {
418 auto& in = r.
inlet(i);
420 if (in.in().type() ==
"Reservoir") {
425 for (
size_t i=0; i<r.
nOutlets(); i++) {
428 if (out.out().type() ==
"Reservoir") {
433 for (
size_t i=0; i<r.
nSurfs(); i++) {
439 if (m_integ ==
nullptr) {
441 "Integrator has not been instantiated. Add one or more reactors first.");
450 m_LHS.assign(m_nv, 1);
451 m_RHS.assign(m_nv, 0);
452 for (
size_t n = 0; n < m_reactors.size(); n++) {
453 m_reactors[n]->applySensitivity(p);
456 if (n == m_reactors.size() - 1) {
461 for (
size_t i =
m_start[n]; i < yEnd; i++) {
462 ydot[i] = m_RHS[i] /
m_LHS[i];
464 m_reactors[n]->resetSensitivity(p);
473 for (
size_t n = 0; n < m_reactors.size(); n++) {
474 m_reactors[n]->applySensitivity(p);
475 m_reactors[n]->evalDae(t, y, ydot, residual);
476 m_reactors[n]->resetSensitivity(p);
483 for (
size_t n = 0; n < m_reactors.size(); n++) {
484 m_reactors[n]->getConstraints(constraints +
m_start[n]);
497 double denom = m_integ->solution(k);
501 return m_integ->sensitivity(k, p) / denom;
508 for (
size_t n = 0; n < m_nv; n++) {
511 double dy = m_atol[n] + fabs(ysave)*m_rtol;
516 eval(t, y, m_ydot.data(), p);
519 for (
size_t m = 0; m < m_nv; m++) {
520 j->
value(m,n) = (m_ydot[m] - ydot[m])/dy;
529 for (
size_t n = 0; n < m_reactors.size(); n++) {
530 m_reactors[n]->updateState(y +
m_start[n]);
539 double* cvode_dky = m_integ->derivative(
m_time, k);
540 for (
size_t j = 0; j < m_nv; j++) {
541 dky[j] = cvode_dky[j];
550 for (
size_t n = 0; n < m_reactors.size(); n++) {
551 m_reactors[n]->setAdvanceLimits(limits +
m_start[n]);
557 bool has_limit =
false;
558 for (
size_t n = 0; n < m_reactors.size(); n++) {
559 has_limit |= m_reactors[n]->hasAdvanceLimits();
566 bool has_limit =
false;
567 for (
size_t n = 0; n < m_reactors.size(); n++) {
568 has_limit |= m_reactors[n]->getAdvanceLimits(limits +
m_start[n]);
575 for (
size_t n = 0; n < m_reactors.size(); n++) {
576 m_reactors[n]->getState(y +
m_start[n]);
582 for (
size_t n = 0; n < m_reactors.size(); n++) {
597 for (
auto r : m_reactors) {
599 return r->name() +
": " + r->componentName(i);
604 throw CanteraError(
"ReactorNet::componentName",
"Index out of bounds");
610 for (
auto r : m_reactors) {
612 return r->upperBound(i);
617 throw CanteraError(
"ReactorNet::upperBound",
"Index {} out of bounds", i0);
623 for (
auto r : m_reactors) {
625 return r->lowerBound(i);
630 throw CanteraError(
"ReactorNet::lowerBound",
"Index {} out of bounds", i0);
635 for (
auto r : m_reactors) {
636 r->resetBadValues(y +
m_start[i++]);
641 const string& name,
double value,
double scale)
644 throw CanteraError(
"ReactorNet::registerSensitivityParameter",
645 "Sensitivity parameters cannot be added after the"
646 "integrator has been initialized.");
657 for (
size_t i = 0; i < m_reactors.size(); i++) {
658 m_reactors[i]->setDerivativeSettings(settings);
665 return m_integ->solverStats();
674 return m_integ->linearSolverType();
684 "Must only be called after ReactorNet is initialized.");
686 m_integ->preconditionerSolve(m_nv, rhs, output);
694 auto precon = m_integ->preconditioner();
698 precon->setGamma(gamma);
700 vector<double> yCopy(m_nv);
704 precon->stateAdjustment(yCopy);
708 for (
size_t i = 0; i < m_reactors.size(); i++) {
709 Eigen::SparseMatrix<double> rJac = m_reactors[i]->jacobian();
710 for (
int k=0; k<rJac.outerSize(); ++k) {
711 for (Eigen::SparseMatrix<double>::InnerIterator it(rJac, k); it; ++it) {
718 precon->updatePreconditioner();
725 "Must only be called after ReactorNet is initialized.");
727 auto precon = m_integ->preconditioner();
728 precon->setGamma(gamma);
729 precon->updatePreconditioner();
734 for (
auto reactor : m_reactors) {
736 throw CanteraError(
"ReactorNet::checkPreconditionerSupported",
737 "Preconditioning is only supported for type *MoleReactor,\n"
738 "Reactor type given: '{}'.",
744SteadyReactorSolver::SteadyReactorSolver(
ReactorNet* net,
double* x0)
747 m_size = m_net->neq();
750 m_initialState.assign(x0, x0 + m_size);
752 m_mask.assign(m_size, 1);
754 for (
size_t i = 0; i < net->nReactors(); i++) {
756 for (
auto& m : R.steadyConstraints()) {
757 m_algebraic.push_back(start + m);
761 for (
auto& n : m_algebraic) {
771 vector<double> xv(x, x +
size());
772 m_net->
eval(0.0, x, r,
nullptr);
773 for (
size_t i = 0; i <
size(); i++) {
791 clock_t t0 = clock();
793 m_work2.resize(
size());
795 for (
size_t j = 0; j <
size(); j++) {
797 double xsave = x0[j];
803 double rdx = 1.0 / (x0[j] - xsave);
806 eval(x0, m_work2.data(), 0.0, 0);
809 for (
size_t i = 0; i <
size(); i++) {
810 double delta = m_work2[i] -
m_work1[i];
812 m_jac->setValue(i, j, delta * rdx);
818 m_jac->updateElapsed(
double(clock() - t0) / CLOCKS_PER_SEC);
819 m_jac->incrementEvals();
826 const double* x =
m_state->data();
827 for (
size_t i = 0; i <
size(); i++) {
828 double ewt = m_net->
rtol()*x[i] + m_net->
atol();
829 double f = step[i] / ewt;
832 return sqrt(sum /
size());
856 const string& message,
int loglevel,
int attempt_counter)
858 if (loglevel >= 6 && !
m_state->empty()) {
860 writelog(
"Current state ({}):\n[", header_suffix);
861 for (
size_t i = 0; i < state.size() - 1; i++) {
866 if (loglevel >= 7 && !
m_xnew.empty()) {
867 writelog(
"Current residual ({}):\n[", header_suffix);
868 for (
size_t i = 0; i <
m_xnew.size() - 1; i++) {
875shared_ptr<ReactorNet>
newReactorNet(vector<shared_ptr<ReactorBase>>& reactors)
877 return make_shared<ReactorNet>(reactors);
Header file for class Cantera::Array2D.
Header file for class ReactorSurface.
Header file for base class WallBase.
A map of string keys to values whose type can vary at runtime.
A class for 2D arrays stored in column-major (Fortran-compatible) form.
double & value(size_t i, size_t j)
Returns a changeable reference to position in the matrix.
Base class for exceptions thrown by Cantera classes.
void setDefaultName(map< string, int > &counts)
Set the default name of a connector. Returns false if it was previously set.
vector< double > m_paramScales
Scaling factors for each sensitivity parameter.
bool suppressErrors() const
Get current state of error suppression.
vector< double > m_sens_params
Values for the problem parameters for which sensitivities are computed This is the array which is per...
An array index is out of range.
Abstract base class for ODE system integrators.
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
virtual void setMaxSteps(int nmax)
Set the maximum number of time-steps the integrator can take before reaching the next output time.
virtual int maxSteps()
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
size_t nWalls()
Return the number of Wall objects connected to this reactor.
WallBase & wall(size_t n)
Return a reference to the n-th Wall connected to this reactor.
void setNetwork(ReactorNet *net)
Set the ReactorNet that this reactor belongs to.
bool setDefaultName(map< string, int > &counts)
Set the default name of a reactor. Returns false if it was previously set.
virtual size_t nSurfs() const
Return the number of surfaces in a reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
size_t nOutlets()
Return the number of outlet FlowDevice objects connected to this reactor.
size_t nInlets()
Return the number of inlet FlowDevice objects connected to this reactor.
ReactorSurface * surface(size_t n)
Return a reference to the n-th ReactorSurface connected to this reactor.
A class representing a network of connected reactors.
void setLinearSolverType(const string &linSolverType="DENSE")
Set the type of linear solver used in the integration.
void preconditionerSetup(double t, double *y, double gamma) override
Evaluate the setup processes for the Jacobian preconditioner.
double step()
Advance the state of all reactors with respect to the independent variable (time or space).
virtual int lastOrder() const
Returns the order used for last solution step of the ODE integrator The function is intended for inte...
void eval(double t, double *y, double *ydot, double *p) override
Evaluate the right-hand-side ODE function.
void initialize()
Initialize the reactor network.
void advance(double t)
Advance the state of all reactors in the independent variable (time or space).
size_t neq() const override
Number of equations.
vector< size_t > m_start
m_start[n] is the starting point in the state vector for reactor n
vector< double > m_LHS
m_LHS is a vector representing the coefficients on the "left hand side" of each governing equation
double m_initial_time
The initial value of the independent variable in the system.
void evalJacobian(double t, double *y, double *ydot, double *p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
void getConstraints(double *constraints) override
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
double m_time
The independent variable in the system.
AnyMap solverStats() const
Get solver stats from integrator.
Reactor & reactor(int n)
Return a reference to the n-th reactor in this network.
map< string, int > m_counts
Map used for default name generation.
double upperBound(size_t i) const
Get the upper bound on the i-th component of the global state vector.
virtual void setMaxSteps(int nmax)
Set the maximum number of internal integration steps the integrator will take before reaching the nex...
string componentName(size_t i) const
Return the name of the i-th component of the global state vector.
void addReactor(Reactor &r)
Add the reactor r to this reactor network.
void getStateDae(double *y, double *ydot) override
Fill in the vectors y and ydot with the current state of the system.
void setInitialTime(double time)
Set the initial value of the independent variable (typically time).
virtual void getDerivative(int k, double *dky)
Return k-th derivative at the current state of the system.
void setMaxErrTestFails(int nmax)
Set the maximum number of error test failures permitted by the CVODES integrator in a single step.
size_t registerSensitivityParameter(const string &name, double value, double scale)
Used by Reactor and Wall objects to register the addition of sensitivity parameters so that the React...
double m_maxstep
Maximum integrator internal timestep. Default of 0.0 means infinity.
double distance()
Current position [m] along the length of the reactor network, for reactors that are solved as a funct...
void setSensitivityTolerances(double rtol, double atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
int maxSteps()
Returns the maximum number of internal integration steps the integrator will take before reaching the...
virtual void setDerivativeSettings(AnyMap &settings)
Set derivative settings of all reactors.
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.
void updateNames(Reactor &r)
Create reproducible names for reactors and walls/connectors.
void updateState(double *y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
void getState(double *y) override
Fill in the vector y with the current state of the system.
void setAdvanceLimits(const double *limits)
Set absolute step size limits during advance.
double rtol()
Relative tolerance.
size_t globalComponentIndex(const string &component, size_t reactor=0)
Return the index corresponding to the component named component in the reactor with index reactor in ...
void solveSteady(int loglevel=0)
Solve directly for the steady-state solution.
bool m_timeIsIndependent
Indicates whether time or space is the independent variable.
double atol()
Absolute integration tolerance.
bool hasAdvanceLimits() const
Check whether ReactorNet object uses advance limits.
void setMaxTimeStep(double maxstep)
Set the maximum integrator step.
double lowerBound(size_t i) const
Get the lower bound on the i-th component of the global state vector.
void evalDae(double t, double *y, double *ydot, double *p, double *residual) override
eval coupling for IDA / DAEs
virtual void checkPreconditionerSupported() const
Check that preconditioning is supported by all reactors in the network.
bool m_integrator_init
True if integrator initialization is current.
void reinitialize()
Reinitialize the integrator.
Integrator & integrator()
Return a reference to the integrator.
bool getAdvanceLimits(double *limits) const
Retrieve absolute step size limits during advance.
Eigen::SparseMatrix< double > steadyJacobian(double rdt=0.0)
Get the Jacobian used by the steady-state solver.
string linearSolverType() const
Problem type of integrator.
void updatePreconditioner(double gamma) override
Update the preconditioner based on already computed jacobian values.
void setPreconditioner(shared_ptr< SystemJacobian > preconditioner)
Set preconditioner used by the linear solver.
void preconditionerSolve(double *rhs, double *output) override
Evaluate the linear system Ax=b where A is the preconditioner.
vector< string > m_paramNames
Names corresponding to each sensitivity parameter.
void resetBadValues(double *y)
Reset physically or mathematically problematic values, such as negative species concentrations.
void setTolerances(double rtol, double atol)
Set the relative and absolute tolerances for the integrator.
virtual void getEstimate(double time, int k, double *yest)
Estimate a future state based on current derivatives.
Class Reactor is a general-purpose class for stirred reactors.
size_t neq()
Number of equations (state variables) for this reactor.
size_t nSensParams() const override
Number of sensitivity parameters associated with this reactor.
string type() const override
String indicating the reactor model implemented.
void initialize(double t0=0.0) override
Initialize the reactor.
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
virtual bool isOde() const
Indicate whether the governing equations for this reactor type are a system of ODEs or DAEs.
virtual bool timeIsIndependent() const
Indicates whether the governing equations for this reactor are functions of time or a spatial variabl...
Adapter class to enable using the SteadyStateSystem solver with ReactorNet.
double weightedNorm(const double *step) const override
Compute the weighted norm of step.
vector< double > m_initialState
Initial value of each state variable.
string componentName(size_t i) const override
Get the name of the i-th component of the state vector.
double upperBound(size_t i) const override
Get the upper bound for global component i in the state vector.
vector< size_t > m_algebraic
Indices of variables that are held constant in the time-stepping mode of the steady-state solver.
void evalJacobian(double *x0) override
Evaluates the Jacobian at x0 using finite differences.
void resetBadValues(double *x) override
Reset values such as negative species concentrations.
void writeDebugInfo(const string &header_suffix, const string &message, int loglevel, int attempt_counter) override
Write solver debugging based on the specified log level.
void eval(double *x, double *r, double rdt=-1.0, int count=1) override
Evaluate the residual function.
void initTimeInteg(double dt, double *x) override
Prepare for time stepping beginning with solution x and timestep dt.
double lowerBound(size_t i) const override
Get the lower bound for global component i in the state vector.
int solve(double *x0, double *x1, int loglevel)
Solve , where is the residual function.
virtual void resize()
Call to set the size of internal data structures after first defining the system or if the problem si...
shared_ptr< SystemJacobian > linearSolver() const
Get the the linear solver being used to hold the Jacobian matrix and solve linear systems as part of ...
vector< double > m_xnew
Work array used to hold the residual or the new solution.
double m_jacobianAbsPerturb
Absolute perturbation of each component in finite difference Jacobian.
size_t size() const
Total solution vector length;.
vector< int > & transientMask()
Access the vector indicating which equations contain a transient term.
double rdt() const
Reciprocal of the time step.
virtual void initTimeInteg(double dt, double *x)
Prepare for time stepping beginning with solution x and timestep dt.
double m_rdt
Reciprocal of time step.
double m_jacobianThreshold
Threshold for ignoring small elements in Jacobian.
shared_ptr< SystemJacobian > m_jac
Jacobian evaluator.
shared_ptr< vector< double > > m_state
Solution vector.
void setMaxTimeStepCount(int nmax)
Set the maximum number of timeteps allowed before successful steady-state solve.
void getState(double *x) const
Get the converged steady-state solution after calling solve().
double m_jacobianRelPerturb
Relative perturbation of each component in finite difference Jacobian.
vector< double > m_work1
Work arrays used during Jacobian evaluation.
void debuglog(const string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Integrator * newIntegrator(const string &itype)
Create new Integrator object.
Namespace for the Cantera kernel.
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
@ BDF_Method
Backward Differentiation.
shared_ptr< SystemJacobian > newSystemJacobian(const string &type)
Create a SystemJacobian object of the specified type.
const double SmallNumber
smallest number to compare to zero.
shared_ptr< ReactorNet > newReactorNet(vector< shared_ptr< ReactorBase > > &reactors)
Create a reactor network containing one or more coupled reactors.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...