1 #ifndef STAN__MATH__ODE__INTEGRATE_ODE_HPP
2 #define STAN__MATH__ODE__INTEGRATE_ODE_HPP
6 #include <boost/numeric/odeint.hpp>
59 template <
typename F,
typename T1,
typename T2>
60 std::vector<std::vector<typename stan::return_type<T1,T2>::type> >
62 const std::vector<T1> y0,
64 const std::vector<double>& ts,
65 const std::vector<T2>& theta,
66 const std::vector<double>& x,
67 const std::vector<int>& x_int,
69 using boost::numeric::odeint::integrate_times;
70 using boost::numeric::odeint::make_dense_output;
71 using boost::numeric::odeint::runge_kutta_dopri5;
74 static_cast<double*>(0));
76 static_cast<double*>(0));
78 static_cast<double*>(0));
80 static_cast<double*>(0));
82 static_cast<double*>(0));
85 static_cast<double*>(0));
87 static_cast<double*>(0));
89 static_cast<double*>(0));
91 static_cast<double*>(0));
93 const double absolute_tolerance = 1
e-6;
94 const double relative_tolerance = 1
e-6;
95 const double step_size = 0.1;
99 coupled_system(f, y0, theta, x, x_int, msgs);
102 std::vector<double> ts_vec(ts.size() + 1);
104 for (
size_t n = 0; n < ts.size(); n++)
107 std::vector<std::vector<double> > y_coupled(ts_vec.size());
111 std::vector<double> initial_coupled_state
112 = coupled_system.initial_state();
114 integrate_times(make_dense_output(absolute_tolerance,
116 runge_kutta_dopri5<std::vector<double>,
121 initial_coupled_state,
122 boost::begin(ts_vec), boost::end(ts_vec),
127 y_coupled.erase(y_coupled.begin());
130 return coupled_system.decouple_states(y_coupled);
bool check_ordered(const char *function, const Eigen::Matrix< T_y, Eigen::Dynamic, 1 > &y, const char *name, T_result *result)
Return true if the specified vector is sorted into increasing order.
bool check_finite(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is finite.
std::vector< std::vector< typename stan::return_type< T1, T2 >::type > > integrate_ode(const F &f, const std::vector< T1 > y0, const 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)
Return the solutions for the specified system of ordinary differential equations given the specified ...
Observer for the coupled states.
double e()
Return the base of the natural logarithm.
Base template class for a coupled ordinary differential equation system, which adds sensitivities to ...
bool check_less(const char *function, const T_y &y, const T_high &high, const char *name, T_result *result)
bool check_nonzero_size(const char *function, const T_y &y, const char *name, T_result *result)
Return true if the specified matrix/vector is of non-zero size.