![]() |
Stan Math Library
5.0.0
Automatic Differentiation
|
Integrator interface for CVODES' adjoint ODE solvers (Adams & BDF methods).
F | Type of ODE right hand side |
T_y0 | Type of scalars for initial state |
T_t0 | Type of initial time |
T_ts | Type of time-points where ODE solution is returned |
T_Args | Types of pass-through parameters |
Definition at line 37 of file cvodes_integrator_adjoint.hpp.
#include <cvodes_integrator_adjoint.hpp>
Classes | |
struct | cvodes_solver |
Since the CVODES solver manages memory with malloc calls, these resources must be freed using a destructor call (which is not being called for the vari class). More... | |
Public Member Functions | |
template<typename FF , require_eigen_col_vector_t< T_y0 > * = nullptr> | |
cvodes_integrator_adjoint_vari (const char *function_name, FF &&f, const T_y0 &y0, const T_t0 &t0, const std::vector< T_ts > &ts, double relative_tolerance_forward, const Eigen::VectorXd &absolute_tolerance_forward, double relative_tolerance_backward, const Eigen::VectorXd &absolute_tolerance_backward, double relative_tolerance_quadrature, double absolute_tolerance_quadrature, long int max_num_steps, long int num_steps_between_checkpoints, int interpolation_polynomial, int solver_forward, int solver_backward, std::ostream *msgs, const T_Args &... args) | |
Construct cvodes_integrator object. | |
std::vector< Eigen::Matrix< T_Return, Eigen::Dynamic, 1 > > | solution () noexcept |
Obtain solution of ODE. | |
void | set_zero_adjoint () final |
No-op for setting adjoints since this class does not own any adjoints. | |
void | chain () final |
Apply the chain rule to this variable based on the variables on which it depends. | |
Static Public Member Functions | |
static void * | operator new (size_t nbytes) noexcept |
Allocate memory from the underlying memory pool. | |
static void | operator delete (void *) noexcept |
Delete a pointer from the underlying memory pool. | |
Private Types | |
using | T_Return = return_type_t< T_y0, T_t0, T_ts, T_Args... > |
using | T_y0_t0 = return_type_t< T_y0, T_t0 > |
Private Member Functions | |
void | store_state (std::size_t n, const Eigen::VectorXd &state, Eigen::Matrix< var, Eigen::Dynamic, 1 > &state_return) |
Overloads which setup the states returned from the forward solve. | |
void | store_state (std::size_t n, const Eigen::VectorXd &state, Eigen::Matrix< double, Eigen::Dynamic, 1 > &state_return) |
template<typename yT , typename... ArgsT> | |
constexpr auto | rhs (double t, const yT &y, const std::tuple< ArgsT... > &args_tuple) const |
Call the ODE RHS with given tuple. | |
int | rhs (double t, const double *y, double *&dy_dt) const |
Calculates the ODE RHS, dy_dt, using the user-supplied functor at the given time t and state y. | |
int | rhs_adj (double t, N_Vector y, N_Vector yB, N_Vector yBdot) const |
int | quad_rhs_adj (double t, N_Vector y, N_Vector yB, N_Vector qBdot) |
int | jacobian_rhs_states (double t, N_Vector y, SUNMatrix J) const |
Calculates the jacobian of the ODE RHS wrt to its states y at the given time-point t and state y. | |
int | jacobian_rhs_adj_states (double t, N_Vector y, SUNMatrix J) const |
Static Private Member Functions | |
static constexpr cvodes_integrator_adjoint_vari * | cast_to_self (void *mem) |
Utility to cast user memory pointer passed in from CVODES to actual typed object pointer. | |
static constexpr int | cv_rhs (realtype t, N_Vector y, N_Vector ydot, void *user_data) |
Implements the function of type CVRhsFn which is the user-defined ODE RHS passed to CVODES. | |
static constexpr int | cv_rhs_adj (realtype t, N_Vector y, N_Vector yB, N_Vector yBdot, void *user_data) |
Implements the function of type CVRhsFnB which is the RHS of the backward ODE system. | |
static constexpr int | cv_quad_rhs_adj (realtype t, N_Vector y, N_Vector yB, N_Vector qBdot, void *user_data) |
Implements the function of type CVQuadRhsFnB which is the RHS of the backward ODE system's quadrature. | |
static constexpr int | cv_jacobian_rhs_states (realtype t, N_Vector y, N_Vector fy, SUNMatrix J, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3) |
Implements the function of type CVDlsJacFn which is the user-defined callback for CVODES to calculate the jacobian of the ode_rhs wrt to the states y. | |
static constexpr int | cv_jacobian_rhs_adj_states (realtype t, N_Vector y, N_Vector yB, N_Vector fyB, SUNMatrix J, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3) |
Implements the CVLsJacFnB function for evaluating the jacobian of the adjoint problem wrt to the backward states. | |
Private Attributes | |
const size_t | num_args_vars_ |
const double | relative_tolerance_forward_ |
const double | relative_tolerance_backward_ |
const double | relative_tolerance_quadrature_ |
const double | absolute_tolerance_quadrature_ |
const long int | max_num_steps_ |
const long int | num_steps_between_checkpoints_ |
const size_t | N_ |
std::ostream * | msgs_ |
vari ** | y_return_varis_ |
vari ** | args_varis_ |
const int | interpolation_polynomial_ |
const int | solver_forward_ |
const int | solver_backward_ |
int | index_backward_ |
bool | backward_is_initialized_ {false} |
cvodes_solver * | solver_ {nullptr} |
Static Private Attributes | |
static constexpr bool | is_var_ts_ {is_var<T_ts>::value} |
static constexpr bool | is_var_t0_ {is_var<T_t0>::value} |
static constexpr bool | is_var_y0_ {is_var<T_y0>::value} |
static constexpr bool | is_var_y0_t0_ {is_var<T_y0_t0>::value} |
static constexpr bool | is_any_var_args_ |
static constexpr bool | is_var_return_ {is_var<T_Return>::value} |
static constexpr bool | is_var_only_ts_ |