Stan Math Library  2.6.3
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerator Friends Macros
integrate_ode.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_ARR_FUNCTOR_INTEGRATE_ODE_HPP
2 #define STAN_MATH_PRIM_ARR_FUNCTOR_INTEGRATE_ODE_HPP
3 
12 #include <boost/numeric/odeint.hpp>
13 #include <ostream>
14 #include <vector>
15 
16 namespace stan {
17 
18  namespace math {
19 
58  template <typename F, typename T1, typename T2>
59  std::vector<std::vector<typename stan::return_type<T1, T2>::type> >
60  integrate_ode(const F& f,
61  const std::vector<T1> y0,
62  const double t0,
63  const std::vector<double>& ts,
64  const std::vector<T2>& theta,
65  const std::vector<double>& x,
66  const std::vector<int>& x_int,
67  std::ostream* msgs) {
68  using boost::numeric::odeint::integrate_times;
69  using boost::numeric::odeint::make_dense_output;
70  using boost::numeric::odeint::runge_kutta_dopri5;
71 
72  stan::math::check_finite("integrate_ode", "initial state", y0);
73  stan::math::check_finite("integrate_ode", "initial time", t0);
74  stan::math::check_finite("integrate_ode", "times", ts);
75  stan::math::check_finite("integrate_ode", "parameter vector", theta);
76  stan::math::check_finite("integrate_ode", "continuous data", x);
77 
78  stan::math::check_nonzero_size("integrate_ode", "times", ts);
79  stan::math::check_nonzero_size("integrate_ode", "initial state", y0);
80  stan::math::check_ordered("integrate_ode", "times", ts);
81  stan::math::check_less("integrate_ode", "initial time", t0, ts[0]);
82 
83  const double absolute_tolerance = 1e-6;
84  const double relative_tolerance = 1e-6;
85  const double step_size = 0.1;
86 
87  // creates basic or coupled system by template specializations
89  coupled_system(f, y0, theta, x, x_int, msgs);
90 
91  // first time in the vector must be time of initial state
92  std::vector<double> ts_vec(ts.size() + 1);
93  ts_vec[0] = t0;
94  for (size_t n = 0; n < ts.size(); n++)
95  ts_vec[n+1] = ts[n];
96 
97  std::vector<std::vector<double> > y_coupled(ts_vec.size());
98  coupled_ode_observer observer(y_coupled);
99 
100  // the coupled system creates the coupled initial state
101  std::vector<double> initial_coupled_state
102  = coupled_system.initial_state();
103 
104  integrate_times(make_dense_output(absolute_tolerance,
105  relative_tolerance,
106  runge_kutta_dopri5<std::vector<double>,
107  double,
108  std::vector<double>,
109  double>() ),
110  coupled_system,
111  initial_coupled_state,
112  boost::begin(ts_vec), boost::end(ts_vec),
113  step_size,
114  observer);
115 
116  // remove the first state corresponding to the initial value
117  y_coupled.erase(y_coupled.begin());
118 
119  // the coupled system also encapsulates the decoupling operation
120  return coupled_system.decouple_states(y_coupled);
121  }
122 
123  }
124 
125 }
126 
127 #endif
bool check_less(const char *function, const char *name, const T_y &y, const T_high &high)
Return true if y is strictly less than high.
Definition: check_less.hpp:81
bool check_nonzero_size(const char *function, const char *name, const T_y &y)
Return true if the specified matrix/vector is of non-zero size.
Observer for the coupled states.
bool check_ordered(const char *function, const char *name, const Eigen::Matrix< T_y, Eigen::Dynamic, 1 > &y)
Return true if the specified vector is sorted into strictly increasing order.
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:95
bool check_finite(const char *function, const char *name, const T_y &y)
Return true if y is finite.
Base template class for a coupled ordinary differential equation system, which adds sensitivities to ...
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 ...

     [ Stan Home Page ] © 2011–2015, Stan Development Team.