27ReactorNet::~ReactorNet()
38 if (reactors.empty()) {
39 throw CanteraError(
"ReactorNet::ReactorNet",
"No reactors in network!");
47 map<Solution*, set<string>> solutions;
51 std::deque<shared_ptr<ReactorBase>> reactorQueue(reactors.begin(), reactors.end());
52 std::set<shared_ptr<ReactorBase>> visited;
54 while (!reactorQueue.empty()) {
55 auto R = reactorQueue.front();
56 reactorQueue.pop_front();
58 if (visited.find(R) != visited.end()) {
63 if (
auto bulk = std::dynamic_pointer_cast<Reactor>(R)) {
64 m_bulkReactors.push_back(bulk);
65 m_reactors.push_back(R);
66 }
else if (
auto surf = std::dynamic_pointer_cast<ReactorSurface>(R)) {
67 m_surfaces.push_back(surf);
68 m_reactors.push_back(R);
69 for (
size_t i = 0; i < surf->nAdjacent(); i++) {
70 reactorQueue.push_back(surf->adjacent(i));
72 }
else if (
auto res = std::dynamic_pointer_cast<Reservoir>(R)) {
73 m_reservoirs.push_back(res);
77 for (
size_t i = 0; i < R->nInlets(); i++) {
78 auto& inlet = R->inlet(i);
79 m_flowDevices.insert(&inlet);
80 reactorQueue.push_back(inlet.in().shared_from_this());
82 for (
size_t i = 0; i < R->nOutlets(); i++) {
83 auto& outlet = R->outlet(i);
84 m_flowDevices.insert(&outlet);
85 reactorQueue.push_back(outlet.out().shared_from_this());
87 for (
size_t i = 0; i < R->nWalls(); i++) {
88 auto& wall = R->wall(i);
89 m_walls.insert(&wall);
90 reactorQueue.push_back(wall.left().shared_from_this());
91 reactorQueue.push_back(wall.right().shared_from_this());
94 auto phase = R->phase();
95 for (
size_t i = 0; i < R->nSurfs(); i++) {
96 if (R->surface(i)->phase()->adjacent(phase->name()) != phase) {
98 "Bulk phase '{}' used by interface '{}' must be the same object\n"
99 "as the contents of the adjacent reactor '{}'.",
100 phase->name(), R->surface(i)->name(), R->name());
102 reactorQueue.push_back(R->surface(i)->shared_from_this());
106 solutions[phase.get()].insert(R->name());
109 for (
auto& R : m_bulkReactors) {
110 if (R->type() ==
"FlowReactor" && m_bulkReactors.size() > 1) {
112 "FlowReactors must be used alone.");
118 size_t nv = R->neq();
121 for (
auto current : m_bulkReactors) {
122 if (current->isOde() != R->isOde()) {
124 "Cannot mix Reactor types using both ODEs and DAEs ({} and {})",
125 current->type(), R->type());
127 if (current->timeIsIndependent() != R->timeIsIndependent()) {
129 "Cannot mix Reactor types using time and space as independent variables"
130 "\n({} and {})", current->type(), R->type());
135 for (
auto surf : m_surfaces) {
136 surf->setOffset(m_nv);
139 solutions[surf->phase().get()].insert(surf->name());
142 for (
auto& [soln, reactors] : solutions) {
143 if (reactors.size() > 1) {
145 for (
auto r : reactors) {
146 if (r != *reactors.begin()) {
149 shared += fmt::format(
"'{}'", r);
151 throw CanteraError(
"ReactorNet::initialize",
"The following reactors /"
152 " reactor surfaces are using the same Solution object: {}. Use"
153 " independent Solution objects or set the 'clone' argument to 'true'"
154 " when creating the Reactor or ReactorSurface objects.", shared);
158 m_ydot.resize(m_nv, 0.0);
159 m_yest.resize(m_nv, 0.0);
160 m_advancelimits.resize(m_nv, -1.0);
167 m_integ->setLinearSolverType(
"DENSE");
197 fill(m_atol.begin(), m_atol.end(), m_atols);
198 m_integ->setTolerances(m_rtol, m_atol);
209 m_integ->setSensitivityTolerances(m_rtolsens, m_atolsens);
216 throw CanteraError(
"ReactorNet::time",
"Time is not the independent variable"
217 " for this reactor network.");
225 throw CanteraError(
"ReactorNet::distance",
"Distance is not the independent"
226 " variable for this reactor network.");
238 if (!m_linearSolverType.empty()) {
239 m_integ->setLinearSolverType(m_linearSolverType);
242 m_integ->setPreconditioner(m_precon);
244 for (
auto&
reactor : m_reactors) {
247 m_integ->initialize(
m_time, *
this);
252 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
262 debuglog(
"Re-initializing reactor network.\n", m_verbose);
263 for (
auto&
reactor : m_reactors) {
266 m_integ->reinitialize(
m_time, *
this);
267 if (m_integ->preconditionerSide() != PreconditionerSide::NO_PRECONDITION) {
278 m_linearSolverType = linSolverType;
284 m_precon = preconditioner;
301 m_integ->integrate(
time);
302 m_time = m_integ->currentTime();
326 m_integ->integrate(
time);
332 m_time = m_integ->currentTime();
346 if (m_advancelimits.size() != m_nv) {
347 m_advancelimits.assign(m_nv, -1.0);
350 auto ycurr = m_integ->solution();
352 double max_ratio = -1.0;
353 double best_limit = 0.0;
354 for (
size_t j = 0; j < m_nv; j++) {
355 double lim = m_advancelimits[j];
357 double delta = std::abs(ycurr[j] -
m_ybase[j]);
358 double ratio = delta / lim;
359 if (ratio > max_ratio) {
368 double y_start =
m_ybase[jmax];
369 double y_end = ycurr[jmax];
370 double delta = y_end - y_start;
371 writelog(
" Advance limit triggered for component {:d} (dt = {:9.4g}):"
372 " y_start = {:11.6g}, y_end = {:11.6g},"
373 " delta = {:11.6g}, limit = {:9.4g}\n",
374 jmax, dt, y_start, y_end, delta, best_limit);
393 vector<double> y(
neq());
397 solver.
solve(loglevel);
405 vector<double> y0(
neq());
406 vector<double> y1(
neq());
413 return std::dynamic_pointer_cast<EigenSparseJacobian>(solver.
linearSolver())->jacobian();
419 auto cvode_dky = m_integ->solution();
420 std::copy(cvode_dky.begin(), cvode_dky.end(), yest.begin());
425 for (
int n = 1; n <= k; n++) {
426 factor *= deltat / n;
427 cvode_dky = m_integ->derivative(
m_time, n);
428 for (
size_t j = 0; j < m_nv; j++) {
429 yest[j] += factor * cvode_dky[j];
437 return m_integ->lastOrder();
455 if (m_advancelimits.size() != m_nv) {
456 m_advancelimits.assign(m_nv, -1.0);
460 double max_ratio = 0.0;
461 for (
size_t i = 0; i < m_nv; i++) {
462 double lim = m_advancelimits[i];
464 double delta = std::abs(y[i] -
m_ybase[i]);
465 double ratio = delta / lim;
466 if (ratio > max_ratio) {
482 for (
size_t i=0; i<r.
nWalls(); i++) {
485 if (w.left().type() ==
"Reservoir") {
488 if (w.right().type() ==
"Reservoir") {
493 for (
size_t i=0; i<r.
nInlets(); i++) {
494 auto& in = r.
inlet(i);
496 if (in.in().type() ==
"Reservoir") {
501 for (
size_t i=0; i<r.
nOutlets(); i++) {
504 if (out.out().type() ==
"Reservoir") {
509 for (
size_t i=0; i<r.
nSurfs(); i++) {
515 if (m_integ ==
nullptr) {
517 "Integrator has not been instantiated. Add one or more reactors first.");
523 span<const double> p)
527 m_LHS.assign(m_nv, 1);
528 m_RHS.assign(m_nv, 0);
529 for (
auto& R : m_reactors) {
530 size_t offset = R->offset();
531 R->applySensitivity(p);
532 R->eval(t, span<double>(
m_LHS.data() +
offset, R->neq()),
533 span<double>(m_RHS.data() +
offset, R->neq()));
535 ydot[i] = m_RHS[i] /
m_LHS[i];
537 R->resetSensitivity(p);
545 m_LHS.assign(m_nv, 1);
546 m_RHS.assign(m_nv, 0);
547 for (
auto& R : m_reactors) {
548 size_t offset = R->offset();
550 span<double>(m_RHS.data() +
offset, R->neq()));
552 residual[i] = m_RHS[i] /
m_LHS[i];
559 span<const double> p, span<double> residual)
563 for (
auto& R : m_reactors) {
564 size_t offset = R->offset();
565 R->applySensitivity(p);
566 R->evalDae(t, y.subspan(
offset, R->neq()),
567 ydot.subspan(
offset, R->neq()),
568 residual.subspan(
offset, R->neq()));
569 R->resetSensitivity(p);
576 for (
auto& R : m_reactors) {
577 R->getConstraints(constraints.subspan(R->offset(), R->neq()));
588 double denom = m_integ->solution(k);
592 return m_integ->sensitivity(k, p) / denom;
596 span<const double> p,
Array2D* j)
600 for (
size_t n = 0; n < m_nv; n++) {
603 double dy = m_atol[n] + fabs(ysave)*m_rtol;
608 eval(t, y, m_ydot, p);
611 for (
size_t m = 0; m < m_nv; m++) {
612 j->
value(m,n) = (m_ydot[m] - ydot[m])/dy;
621 for (
auto& R : m_reactors) {
622 R->updateState(y.subspan(R->offset(), R->neq()));
629 auto cvode_dky = m_integ->derivative(
m_time, k);
630 std::copy(cvode_dky.begin(), cvode_dky.end(), dky.begin());
636 for (
auto& R : m_bulkReactors) {
637 R->setAdvanceLimits(limits.subspan(R->offset(), R->neq()));
643 bool has_limit =
false;
644 for (
size_t n = 0; n < m_bulkReactors.size(); n++) {
645 has_limit |= m_bulkReactors[n]->hasAdvanceLimits();
652 bool has_limit =
false;
653 for (
auto& R : m_bulkReactors) {
654 has_limit |= R->getAdvanceLimits(limits.subspan(R->offset(), R->neq()));
661 for (
auto& R : m_reactors) {
662 R->getState(y.subspan(R->offset(), R->neq()));
672 for (
size_t n = m_reactors.size(); n != 0 ; n--) {
673 auto& R = m_reactors[n-1];
674 R->getStateDae(y.subspan(R->offset(), R->neq()),
675 ydot.subspan(R->offset(), R->neq()));
682 return m_reactors[
reactor]->offset()
683 + m_reactors[
reactor]->componentIndex(component);
690 for (
auto r : m_reactors) {
692 return r->name() +
": " + r->componentName(i);
698 throw IndexError(
"ReactorNet::componentName",
"component", i0, iTot);
705 for (
auto r : m_reactors) {
707 return r->upperBound(i);
713 throw IndexError(
"ReactorNet::upperBound",
"upperBound", i0, iTot);
720 for (
auto r : m_reactors) {
722 return r->lowerBound(i);
728 throw IndexError(
"ReactorNet::lowerBound",
"lowerBound", i0, iTot);
732 for (
auto& R : m_reactors) {
733 R->resetBadValues(y.subspan(R->offset(), R->neq()));
738 const string& name,
double value,
double scale)
741 throw CanteraError(
"ReactorNet::registerSensitivityParameter",
742 "Sensitivity parameters cannot be added after the"
743 "integrator has been initialized.");
754 for (
auto& R : m_bulkReactors) {
755 R->setDerivativeSettings(settings);
762 return m_integ->solverStats();
771 return m_integ->linearSolverType();
781 "Must only be called after ReactorNet is initialized.");
783 m_integ->preconditionerSolve(rhs, output);
791 auto precon = std::dynamic_pointer_cast<EigenSparseJacobian>(
792 m_integ->preconditioner());
796 precon->setGamma(gamma);
798 vector<double> yCopy(m_nv);
802 precon->stateAdjustment(yCopy);
806 vector<Eigen::Triplet<double>> trips;
807 for (
auto& R : m_reactors) {
808 R->getJacobianElements(trips);
810 precon->setFromTriplets(trips);
812 precon->updatePreconditioner();
819 "Must only be called after ReactorNet is initialized.");
821 auto precon = m_integ->preconditioner();
822 precon->setGamma(gamma);
823 precon->updatePreconditioner();
828 for (
auto reactor : m_bulkReactors) {
830 throw CanteraError(
"ReactorNet::checkPreconditionerSupported",
831 "Preconditioning is only supported for type *MoleReactor,\n"
832 "Reactor type given: '{}'.",
838SteadyReactorSolver::SteadyReactorSolver(
ReactorNet* net, span<const double> x0)
841 m_size = m_net->neq();
844 m_initialState.assign(x0.begin(), x0.end());
846 setJacobianPerturbation(m_jacobianRelPerturb, 1000 * m_net->atol(),
847 m_jacobianThreshold);
848 m_mask.assign(m_size, 1);
850 for (
size_t i = 0; i < net->nReactors(); i++) {
853 auto algebraic = R.initializeSteady();
854 for (
auto& m : algebraic) {
855 m_mask[start + m] = 0;
862 double rdt,
int count)
868 for (
size_t i = 0; i <
size(); i++) {
882 clock_t t0 = clock();
884 m_work2.resize(
size());
886 vector<double> perturbed(x0.begin(), x0.end());
887 for (
size_t j = 0; j <
size(); j++) {
889 double xsave = x0[j];
894 perturbed[j] = xsave + dx;
895 double rdx = 1.0 / (perturbed[j] - xsave);
898 eval(perturbed, m_work2, 0.0, 0);
901 for (
size_t i = 0; i <
size(); i++) {
902 double delta = m_work2[i] -
m_work1[i];
904 m_jac->setValue(i, j, delta * rdx);
907 perturbed[j] = xsave;
912 m_jac->updateElapsed(
double(clock() - t0) / CLOCKS_PER_SEC);
913 m_jac->incrementEvals();
920 const double* x =
m_state->data();
921 for (
size_t i = 0; i <
size(); i++) {
922 double ewt = m_net->
rtol()*fabs(x[i]) + m_net->
atol();
923 double f = step[i] / ewt;
926 return sqrt(sum /
size());
950 const string& message,
int loglevel,
int attempt_counter)
952 if (loglevel >= 6 && !
m_state->empty()) {
954 writelog(
"Current state ({}):\n[", header_suffix);
955 for (
size_t i = 0; i < state.size() - 1; i++) {
960 if (loglevel >= 7 && !
m_xnew.empty()) {
961 writelog(
"Current residual ({}):\n[", header_suffix);
962 for (
size_t i = 0; i <
m_xnew.size() - 1; i++) {
971 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.
Base class for reactor objects.
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.
virtual string type() const
String indicating the reactor model implemented.
bool setDefaultName(map< string, int > &counts)
Set the default name of a reactor. Returns false if it was previously set.
virtual void updateState(span< const double > y)
Set the state of the reactor to correspond to the state vector y.
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.
virtual bool preconditionerSupported() const
Return a false if preconditioning is not supported or true otherwise.
virtual void updateConnected(bool updatePressure)
Update state information needed by connected reactors, flow devices, and walls.
A class representing a network of connected reactors.
virtual void getDerivative(int k, span< double > dky)
Return k-th derivative at the current state of the system.
void setLinearSolverType(const string &linSolverType="DENSE")
Set the type of linear solver used in the integration.
bool getAdvanceLimits(span< double > limits) const
Retrieve absolute step size limits during advance.
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...
double m_ybase_time
Base time corresponding to m_ybase.
void initialize()
Initialize the reactor network.
void evalDae(double t, span< const double > y, span< const double > ydot, span< const double > p, span< double > residual) override
eval coupling for IDA / DAEs
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.
ReactorBase & reactor(int n)
Return a reference to the n-th reactor in this network.
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.
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
void updateNames(ReactorBase &r)
Create reproducible names for reactors and walls/connectors.
double m_time
The independent variable in the system.
AnyMap solverStats() const
Get solver stats from integrator.
bool m_limit_check_active
Indicates whether the advance-limit root check is active for the current call to advance(t,...
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...
bool m_integratorInitialized
True if the integrator has been initialized at least once.
string componentName(size_t i) const
Return the name of the i-th component of the global state vector.
virtual void getEstimate(double time, int k, span< double > yest)
Estimate a future state based on current derivatives.
size_t nRootFunctions() const override
Root finding is enabled only while enforcing advance limits.
void evalRootFunctions(double t, span< const double > y, span< double > gout) override
Evaluate the advance-limit root function used to stop integration once a limit is met.
void setInitialTime(double time)
Set the initial value of the independent variable (typically time).
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...
void getConstraints(span< double > constraints) override
Given a vector of length neq(), mark which variables should be considered algebraic constraints.
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 evalSteady(span< const double > y, span< double > residual)
Evaluate the governing equations adapted for the steady-state solver.
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.
ReactorNet(shared_ptr< ReactorBase > reactor)
Create reactor network containing single reactor.
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 evalJacobian(double t, span< double > y, span< double > ydot, span< const double > p, Array2D *j)
Evaluate the Jacobian matrix for the reactor network.
void updateState(span< const double > y)
Update the state of all the reactors in the network to correspond to the values in the solution vecto...
void preconditionerSolve(span< const double > rhs, span< double > output) override
Evaluate the linear system Ax=b where A is the preconditioner.
bool m_needIntegratorInit
True if integrator needs to be (re)initialized.
void preconditionerSetup(double t, span< const double > y, double gamma) override
Evaluate the setup processes for the Jacobian preconditioner.
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 resetBadValues(span< double > y)
Reset physically or mathematically problematic values, such as negative species concentrations.
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.
vector< double > m_ybase
Base state used for evaluating advance limits during a single advance() call when root-finding is ena...
void setMaxTimeStep(double maxstep)
Set the maximum integrator step.
void setAdvanceLimits(span< const double > limits)
Set absolute step size limits during advance.
double lowerBound(size_t i) const
Get the lower bound on the i-th component of the global state vector.
virtual void checkPreconditionerSupported() const
Check that preconditioning is supported by all reactors in the network.
void reinitialize()
Reinitialize the integrator.
Integrator & integrator()
Return a reference to the integrator.
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 getState(span< double > y) override
Fill in the vector y with the current state of the system.
vector< string > m_paramNames
Names corresponding to each sensitivity parameter.
void getStateDae(span< double > y, span< double > ydot) override
Fill in the vectors y and ydot with the current state of the system.
void setTolerances(double rtol, double atol)
Set the relative and absolute tolerances for the integrator.
void eval(double t, span< const double > y, span< double > ydot, span< const double > p) override
Evaluate the right-hand-side ODE function.
Adapter class to enable using the SteadyStateSystem solver with ReactorNet.
double weightedNorm(span< 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.
void initTimeInteg(double dt, span< const double > x) override
Prepare for time stepping beginning with solution x and timestep dt.
double upperBound(size_t i) const override
Get the upper bound for global component i in the state vector.
void evalJacobian(span< const double > x0) override
Evaluates the Jacobian at x0 using finite differences.
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 resetBadValues(span< double > x) override
Reset values such as negative species concentrations.
void eval(span< const double > x, span< double > r, double rdt=-1.0, int count=1) override
Evaluate the residual function.
double lowerBound(size_t i) const override
Get the lower bound for global component i in the state vector.
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.
void getState(span< double > x) const
Get the converged steady-state solution after calling solve().
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, span< const 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.
vector< int > m_mask
Transient mask.
void setMaxTimeStepCount(int nmax)
Set the maximum number of timeteps allowed before successful steady-state solve.
void solve(int loglevel=0)
Solve the steady-state problem, taking internal timesteps as necessary until the Newton solver can co...
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.
const size_t npos
index returned by functions to indicate "no position"
void checkFinite(const double tmp)
Check to see that a number is finite (not NaN, +Inf or -Inf)
shared_ptr< ReactorNet > newReactorNet(span< shared_ptr< ReactorBase > > reactors)
Create a reactor network containing one or more coupled reactors.
@ 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.
offset
Offsets of solution components in the 1D solution array.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...