Class CVodeSolver
Defined in File solver_cvodes.h
Inheritance Relationships
Base Type
public amici::Solver(Class Solver)
Class Documentation
-
class CVodeSolver : public amici::Solver
The CVodeSolver class is a wrapper around the SUNDIALS CVODES solver.
Public Functions
-
~CVodeSolver() override = default
-
inline virtual std::string get_class_name() const override
Get the name of this class.
- Returns:
Class name.
-
virtual void reinit(realtype t0, AmiVector const &yy0, AmiVector const &yp0) const override
Reinitializes the states in the solver after an event occurrence.
- Parameters:
t0 – reinitialization timepoint
yy0 – initial state variables
yp0 – initial derivative state variables (DAE only)
-
virtual void sens_reinit(AmiVectorArray const &yyS0, AmiVectorArray const &ypS0) const override
Reinitializes the state sensitivities in the solver after an event occurrence.
- Parameters:
yyS0 – new state sensitivity
ypS0 – new derivative state sensitivities (DAE only)
-
virtual void sens_toggle_off() const override
Switches off computation of state sensitivities without deallocating the memory for sensitivities.
-
virtual void reinit_b(int which, realtype tB0, AmiVector const &yyB0, AmiVector const &ypB0) const override
Reinitializes the adjoint states after an event occurrence.
- Parameters:
which – identifier of the backwards problem
tB0 – reinitialization timepoint
yyB0 – new adjoint state
ypB0 – new adjoint derivative state
-
virtual void reinit_quad_b(int which, AmiVector const &yQB0) const override
Reinitialize the adjoint states after an event occurrence.
- Parameters:
which – identifier of the backwards problem
yQB0 – new adjoint quadrature state
-
virtual int solve(realtype tout, int itask) const override
Solves the forward problem until a predefined timepoint.
- Parameters:
tout – timepoint until which simulation should be performed
itask – task identifier, can be CV_NORMAL or CV_ONE_STEP
- Returns:
status flag indicating success of execution
-
virtual int solve_f(realtype tout, int itask, int *ncheckPtr) const override
Solves the forward problem until a predefined timepoint (adjoint only)
- Parameters:
tout – timepoint until which simulation should be performed
itask – task identifier, can be CV_NORMAL or CV_ONE_STEP
ncheckPtr – pointer to a number that counts the internal checkpoints
- Returns:
status flag indicating success of execution
-
virtual void solve_b(realtype tBout, int itaskB) const override
Solves the backward problem until a predefined timepoint (adjoint only)
- Parameters:
tBout – timepoint until which simulation should be performed
itaskB – task identifier, can be CV_NORMAL or CV_ONE_STEP
-
virtual void get_dky(realtype t, int k) const override
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
-
virtual void get_sens_dky(realtype t, int k) const override
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
-
virtual void get_quad_dky_b(realtype t, int k, int which) const override
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
which – index of backward problem
-
virtual void get_dky_b(realtype t, int k, int which) const override
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
which – index of backward problem
-
virtual void get_root_info(int *rootsfound) const override
get_root_info extracts information which event occurred
- Parameters:
rootsfound – array with flags indicating whether the respective event occurred
-
virtual void set_stop_time(realtype tstop) const override
Sets a timepoint at which the simulation will be stopped.
- Parameters:
tstop – timepoint until which simulation should be performed
-
virtual void turn_off_root_finding() const override
Disable rootfinding.
-
virtual Model const *get_model() const override
Accessor function to the model stored in the user data
- Returns:
user data model
-
virtual void set_linear_solver() const override
Sets the linear solver for the forward problem.
-
virtual void set_linear_solver_b(int which) const override
Sets the linear solver for the backward problem.
- Parameters:
which – index of the backward problem
-
virtual void set_non_linear_solver() const override
Set the non-linear solver for the forward problem.
-
virtual void set_non_linear_solver_sens() const override
Set the non-linear solver for sensitivities.
-
virtual void set_non_linear_solver_b(int which) const override
Set the non-linear solver for the backward problem.
- Parameters:
which – index of the backward problem
-
Solver() = default
Default constructor.
Protected Functions
-
virtual void calc_ic(realtype tout1) const override
Calculates consistent initial conditions, assumes initial states to be correct (DAE only)
- Parameters:
tout1 – next timepoint to be computed (sets timescale)
-
virtual void calc_ic_b(int which, realtype tout1) const override
Calculates consistent initial conditions for the backwards problem, assumes initial states to be correct (DAE only)
- Parameters:
which – identifier of the backwards problem
tout1 – next timepoint to be computed (sets timescale)
-
virtual void get_b(int which) const override
extracts the adjoint state at the current timepoint from solver memory and writes it to the xB member variable
- Parameters:
which – index of the backwards problem
-
virtual void get_sens() const override
extracts the state sensitivity at the current timepoint from solver memory and writes it to the sx member variable
-
virtual void get_quad_b(int which) const override
extracts the adjoint quadrature state at the current timepoint from solver memory and writes it to the xQB member variable
- Parameters:
which – index of the backwards problem
-
virtual void get_quad(realtype &t) const override
extracts the quadrature at the current timepoint from solver memory and writes it to the xQ member variable
- Parameters:
t – timepoint for quadrature extraction
-
virtual void get_quad_dky(realtype t, int k) const override
interpolates the (derivative of the) solution at the requested timepoint
- Parameters:
t – timepoint
k – derivative order
-
virtual void reinit_post_process_f(realtype tnext) const override
reinit_post_process_f postprocessing of the solver memory after a discontinuity in the forward problem
- Parameters:
tnext – next timepoint (defines integration direction)
-
virtual void reinit_post_process_b(realtype tnext) const override
reinit_post_process_b postprocessing of the solver memory after a discontinuity in the backward problem
- Parameters:
tnext – next timepoint (defines integration direction)
-
void reInit_post_process(void *cv_mem, realtype *t, AmiVector *yout, realtype tout) const
Post-processing of the solver memory after a discontinuity.
- Parameters:
cv_mem – pointer to CVODES solver memory object
t – pointer to integration time
yout – new state vector
tout – anticipated next integration timepoint.
-
virtual void allocate_solver() const override
Create specifies solver method and initializes solver memory for the forward problem.
-
virtual void set_ss_tolerances(double rtol, double atol) const override
sets scalar relative and absolute tolerances for the forward problem
- Parameters:
rtol – relative tolerances
atol – absolute tolerances
-
virtual void set_sens_ss_tolerances(double rtol, double const *atol) const override
activates sets scalar relative and absolute tolerances for the sensitivity variables
- Parameters:
rtol – relative tolerances
atol – array of absolute tolerances for every sensitivity variable
-
virtual void set_sens_err_con(bool error_corr) const override
SetSensErrCon specifies whether error control is also enforced for sensitivities for the forward problem
- Parameters:
error_corr – activation flag
-
virtual void set_quad_err_con_b(int which, bool flag) const override
Specifies whether error control is also enforced for the backward quadrature problem.
- Parameters:
which – identifier of the backwards problem
flag – activation flag
-
virtual void set_quad_err_con(bool flag) const override
Specifies whether error control is also enforced for the forward quadrature problem.
- Parameters:
flag – activation flag
-
virtual void set_user_data() const override
Attaches the user data to the forward problem.
-
virtual void set_user_data_b(int which) const override
attaches the user data to the backward problem
- Parameters:
which – identifier of the backwards problem
-
virtual void set_max_num_steps(long int mxsteps) const override
specifies the maximum number of steps for the forward problem
Note
in contrast to the SUNDIALS method, this sets the overall maximum, not the maximum between output times.
- Parameters:
mxsteps – number of steps
-
virtual void set_stab_lim_det(int stldet) const override
activates stability limit detection for the forward problem
- Parameters:
stldet – flag for stability limit detection (TRUE or FALSE)
-
virtual void set_stab_lim_det_b(int which, int stldet) const override
activates stability limit detection for the backward problem
- Parameters:
which – identifier of the backwards problem
stldet – flag for stability limit detection (TRUE or FALSE)
-
virtual void set_id(Model const *model) const override
specify algebraic/differential components (DAE only)
- Parameters:
model – model specification
-
virtual void set_suppress_alg(bool flag) const override
deactivates error control for algebraic components (DAE only)
- Parameters:
flag – deactivation flag
-
void reset_state(void *cv_mem, const_N_Vector y0) const
reset_state reset the CVODES solver to restart integration after a rhs discontinuity.
- Parameters:
cv_mem – pointer to CVODES solver memory object
y0 – new state vector
-
virtual void set_sens_params(realtype const *p, realtype const *pbar, int const *plist) const override
specifies the scaling and indexes for sensitivity computation
- Parameters:
p – parameters
pbar – parameter scaling constants
plist – parameter index list
-
virtual void adj_init() const override
initializes the adjoint problem
-
virtual void quad_init(AmiVector const &xQ0) const override
initializes the quadratures
- Parameters:
xQ0 – vector with initial values for xQ
-
virtual void allocate_solver_b(int *which) const override
Specifies solver method and initializes solver memory for the backward problem.
- Parameters:
which – identifier of the backwards problem
-
virtual void set_ss_tolerances_b(int which, realtype relTolB, realtype absTolB) const override
sets relative and absolute tolerances for the backward problem
- Parameters:
which – identifier of the backwards problem
relTolB – relative tolerances
absTolB – absolute tolerances
-
virtual void quad_ss_tolerances_b(int which, realtype reltolQB, realtype abstolQB) const override
sets relative and absolute tolerances for the quadrature backward problem
- Parameters:
which – identifier of the backwards problem
reltolQB – relative tolerances
abstolQB – absolute tolerances
-
virtual void quad_ss_tolerances(realtype reltolQ, realtype abstolQ) const override
sets relative and absolute tolerances for the quadrature problem
- Parameters:
reltolQB – relative tolerances
abstolQB – absolute tolerances
-
virtual void set_max_num_steps_b(int which, long int mxstepsB) const override
specifies the maximum number of steps for the forward problem
Note
in contrast to the SUNDIALS method, this sets the overall maximum, not the maximum between output times.
- Parameters:
which – identifier of the backwards problem
mxstepsB – number of steps
-
virtual void diag() const override
attaches a diagonal linear solver to the forward problem
-
virtual void diag_b(int which) const override
attaches a diagonal linear solver to the backward problem
- Parameters:
which – identifier of the backwards problem
-
virtual void get_num_steps(void const *ami_mem, long int *numsteps) const override
reports the number of solver steps
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
numsteps – output array
-
virtual void get_num_rhs_evals(void const *ami_mem, long int *numrhsevals) const override
reports the number of right hand evaluations
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
numrhsevals – output array
-
virtual void get_num_err_test_fails(void const *ami_mem, long int *numerrtestfails) const override
reports the number of local error test failures
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
numerrtestfails – output array
-
virtual void get_num_non_lin_solv_conv_fails(void const *ami_mem, long int *numnonlinsolvconvfails) const override
reports the number of nonlinear convergence failures
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
numnonlinsolvconvfails – output array
-
virtual void get_last_order(void const *ami_mem, int *order) const override
Reports the order of the integration method during the last internal step.
- Parameters:
ami_mem – pointer to the solver memory instance (can be from forward or backward problem)
order – output array
-
virtual void *get_adj_b_mem(void *ami_mem, int which) const override
Retrieves the solver memory instance for the backward problem.
- Parameters:
which – identifier of the backwards problem
ami_mem – pointer to the forward solver memory instance
- Returns:
A (void *) pointer to the CVODES memory allocated for the backward problem.
-
virtual void init(realtype t0, AmiVector const &x0, AmiVector const &dx0) const override
Initializes the states at the specified initial timepoint.
- Parameters:
t0 – initial timepoint
x0 – initial states
dx0 – initial derivative states
-
virtual void init_steady_state(realtype t0, AmiVector const &x0, AmiVector const &dx0) const override
Initializes the states at the specified initial timepoint.
- Parameters:
t0 – initial timepoint
x0 – initial states
dx0 – initial derivative states
-
virtual void sens_init_1(AmiVectorArray const &sx0, AmiVectorArray const &sdx0) const override
Initializes the forward sensitivities.
- Parameters:
sx0 – initial states sensitivities
sdx0 – initial derivative states sensitivities
-
virtual void b_init(int which, realtype tf, AmiVector const &xB0, AmiVector const &dxB0) const override
Initialize the adjoint states at the specified final timepoint.
- Parameters:
which – identifier of the backwards problem
tf – final timepoint
xB0 – initial adjoint state
dxB0 – initial adjoint derivative state
-
virtual void qb_init(int which, AmiVector const &xQB0) const override
Initialize the quadrature states at the specified final timepoint.
- Parameters:
which – identifier of the backwards problem
xQB0 – initial adjoint quadrature state
-
virtual void root_init(int ne) const override
Initializes the root-finding for events.
- Parameters:
ne – number of different events
-
virtual void set_dense_jac_fn() const override
Set the dense Jacobian function.
-
virtual void set_sparse_jac_fn() const override
sets the sparse Jacobian function
-
virtual void set_band_jac_fn() const override
sets the banded Jacobian function
-
virtual void set_jac_times_vec_fn() const override
sets the Jacobian vector multiplication function
-
virtual void set_dense_jac_fn_b(int which) const override
sets the dense Jacobian function
- Parameters:
which – identifier of the backwards problem
-
virtual void set_sparse_jac_fn_b(int which) const override
sets the sparse Jacobian function
- Parameters:
which – identifier of the backwards problem
-
virtual void set_band_jac_fn_b(int which) const override
sets the banded Jacobian function
- Parameters:
which – identifier of the backwards problem
-
virtual void set_jac_times_vec_fn_b(int which) const override
sets the Jacobian vector multiplication function
- Parameters:
which – identifier of the backwards problem
-
virtual void set_sparse_jac_fn_ss() const override
sets the sparse Jacobian function for backward steady state case
-
virtual void apply_max_nonlin_iters() const override
Apply the maximum number of nonlinear solver iterations permitted per step.
-
virtual void apply_max_conv_fails() const override
Apply the maximum number of nonlinear solver convergence failures permitted per step.
-
virtual void apply_constraints() const override
Apply the constraints to the solver.
-
virtual void apply_max_step_size() const override
Apply the allowed maximum stepsize to the solver.
Friends
-
template<class Archive>
friend void serialize(Archive &ar, CVodeSolver &s, unsigned int) Serialize amici::CVodeSolver to boost archive.
- Parameters:
ar – Archive
s – Solver instance to serialize
-
friend bool operator==(CVodeSolver const &a, CVodeSolver const &b)
Equality operator.
- Parameters:
a –
b –
- Returns:
Whether a and b are equal
-
~CVodeSolver() override = default