6 #ifndef CT_REACTORNET_H
7 #define CT_REACTORNET_H
18 class PreconditionerBase;
132 double advance(
double t,
bool applylimit);
145 return *m_reactors[n];
205 double* ydot,
double* p,
Array2D* j);
208 size_t neq()
const override {
212 size_t nReactors()
const {
213 return m_reactors.size();
216 void eval(
double t,
double* y,
double* ydot,
double* p)
override;
219 void evalDae(
double t,
double* y,
double* ydot,
double* p,
220 double* residual)
override;
223 void getStateDae(
double* y,
double* ydot)
override;
318 vector<Reactor*> m_reactors;
319 unique_ptr<Integrator> m_integ;
335 vector<double> m_atol;
336 double m_rtol = 1.0e-9;
337 double m_rtolsens = 1.0e-4;
338 double m_atols = 1.0e-15;
339 double m_atolsens = 1.0e-6;
340 shared_ptr<PreconditionerBase> m_precon;
341 string m_linearSolverType;
346 bool m_verbose =
false;
354 vector<double> m_ydot;
355 vector<double> m_yest;
356 vector<double> m_advancelimits;
360 vector<double> m_RHS;
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.
Virtual base class for ODE/DAE right-hand-side function evaluators.
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...
Abstract base class for ODE system integrators.
A class representing a network of connected reactors.
void setPreconditioner(shared_ptr< PreconditionerBase > preconditioner)
Set preconditioner used by the linear solver.
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.
size_t nparams() const override
Number of sensitivity parameters.
void initialize()
Initialize the reactor network.
void advance(double t)
Advance the state of all reactors in the independent variable (time or space).
const string & sensitivityParameterName(size_t p) const
The name of the p-th sensitivity parameter added to this ReactorNet.
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.
double getInitialTime() const
Get the initial value of the independent variable (typically time).
void setNeedsReinit()
Called to trigger integrator reinitialization before further integration.
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.
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 maxTimeStep() const
Get the maximum integrator step.
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...
double atolSensitivity() const
Absolute sensitivity tolerance.
void setSensitivityTolerances(double rtol, double atol)
Set the relative and absolute tolerances for integrating the sensitivity equations.
Reactor & reactor(int n)
Return a reference to the n-th reactor in this network.
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 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.
bool verbose() const
Returns true if verbose logging output is enabled.
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 ...
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.
double sensitivity(const string &component, size_t p, int reactor=0)
Return the sensitivity of the component named component with respect to the p-th sensitivity paramete...
void setMaxTimeStep(double maxstep)
Set the maximum integrator step.
double rtolSensitivity() const
Relative sensitivity tolerance.
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.
void setVerbose(bool v=true)
Enable or disable verbose logging while setting up and integrating the reactor network.
string linearSolverType() const
Problem type of integrator.
void updatePreconditioner(double gamma) override
Update the preconditioner based on already computed jacobian values.
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 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.
void scale(InputIter begin, InputIter end, OutputIter out, S scale_factor)
Multiply elements of an array by a scale factor.
Namespace for the Cantera kernel.