1 #ifndef STAN_MATH_PRIM_ARR_FUNCTOR_INTEGRATE_ODE_RK45_HPP 2 #define STAN_MATH_PRIM_ARR_FUNCTOR_INTEGRATE_ODE_RK45_HPP 13 #include <boost/numeric/odeint.hpp> 64 template <
typename F,
typename T1,
typename T2>
65 std::vector<std::vector<typename stan::return_type<T1, T2>::type> >
67 const std::vector<T1> y0,
69 const std::vector<double>& ts,
70 const std::vector<T2>& theta,
71 const std::vector<double>& x,
72 const std::vector<int>& x_int,
73 std::ostream* msgs = 0,
74 double relative_tolerance = 1
e-6,
75 double absolute_tolerance = 1
e-6,
76 int max_num_steps = 1E6) {
77 using boost::numeric::odeint::integrate_times;
78 using boost::numeric::odeint::make_dense_output;
79 using boost::numeric::odeint::runge_kutta_dopri5;
80 using boost::numeric::odeint::max_step_checker;
85 check_finite(
"integrate_ode_rk45",
"parameter vector", theta);
86 check_finite(
"integrate_ode_rk45",
"continuous data", x);
91 check_less(
"integrate_ode_rk45",
"initial time", t0, ts[0]);
93 if (relative_tolerance <= 0)
95 "relative_tolerance,", relative_tolerance,
96 "",
", must be greater than 0");
97 if (absolute_tolerance <= 0)
99 "absolute_tolerance,", absolute_tolerance,
100 "",
", must be greater than 0");
101 if (max_num_steps <= 0)
103 "max_num_steps,", max_num_steps,
104 "",
", must be greater than 0");
108 coupled_system(f, y0, theta, x, x_int, msgs);
111 std::vector<double> ts_vec(ts.size() + 1);
113 for (
size_t n = 0; n < ts.size(); n++)
116 std::vector<std::vector<double> > y_coupled(ts_vec.size());
120 std::vector<double> initial_coupled_state
121 = coupled_system.initial_state();
123 const double step_size = 0.1;
124 integrate_times(make_dense_output(absolute_tolerance,
126 runge_kutta_dopri5<std::vector<double>,
131 initial_coupled_state,
132 boost::begin(ts_vec), boost::end(ts_vec),
135 max_step_checker(max_num_steps));
138 y_coupled.erase(y_coupled.begin());
141 return coupled_system.decouple_states(y_coupled);
void check_finite(const char *function, const char *name, const T_y &y)
Check if y is finite.
void check_nonzero_size(const char *function, const char *name, const T_y &y)
Check if the specified matrix/vector is of non-zero size.
void check_ordered(const char *function, const char *name, const std::vector< T_y > &y)
Check if the specified vector is sorted into strictly increasing order.
std::vector< std::vector< typename stan::return_type< T1, T2 >::type > > integrate_ode_rk45(const F &f, const std::vector< T1 > y0, double t0, const std::vector< double > &ts, const std::vector< T2 > &theta, const std::vector< double > &x, const std::vector< int > &x_int, std::ostream *msgs=0, double relative_tolerance=1e-6, double absolute_tolerance=1e-6, int max_num_steps=1E6)
Return the solutions for the specified system of ordinary differential equations given the specified ...
Observer for the coupled states.
void invalid_argument(const char *function, const char *name, const T &y, const char *msg1, const char *msg2)
Throw an invalid_argument exception with a consistently formatted message.
double e()
Return the base of the natural logarithm.
Base template class for a coupled ordinary differential equation system, which adds sensitivities to ...
void check_less(const char *function, const char *name, const T_y &y, const T_high &high)
Check if y is strictly less than high.