Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
integrate_ode.hpp
Go to the documentation of this file.
1 #ifndef STAN__MATH__ODE__INTEGRATE_ODE_HPP
2 #define STAN__MATH__ODE__INTEGRATE_ODE_HPP
3 
4 #include <ostream>
5 #include <vector>
6 #include <boost/numeric/odeint.hpp>
7 #include <stan/meta/traits.hpp>
13 
16 
17 namespace stan {
18 
19  namespace math {
20 
59  template <typename F, typename T1, typename T2>
60  std::vector<std::vector<typename stan::return_type<T1,T2>::type> >
61  integrate_ode(const F& f,
62  const std::vector<T1> y0,
63  const double t0,
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,
68  std::ostream* msgs) {
69  using boost::numeric::odeint::integrate_times;
70  using boost::numeric::odeint::make_dense_output;
71  using boost::numeric::odeint::runge_kutta_dopri5;
72 
73  stan::math::check_finite("integrate_ode(%1%)", y0, "initial state",
74  static_cast<double*>(0));
75  stan::math::check_finite("integrate_ode(%1%)", t0, "initial time",
76  static_cast<double*>(0));
77  stan::math::check_finite("integrate_ode(%1%)", ts, "times",
78  static_cast<double*>(0));
79  stan::math::check_finite("integrate_ode(%1%)", theta, "parameter vector",
80  static_cast<double*>(0));
81  stan::math::check_finite("integrate_ode(%1%)", x, "continuous data",
82  static_cast<double*>(0));
83 
84  stan::math::check_nonzero_size("integrate_ode(%1%)", ts, "times",
85  static_cast<double*>(0));
86  stan::math::check_nonzero_size("integrate_ode(%1%)", y0, "initial state",
87  static_cast<double*>(0));
88  stan::math::check_ordered("integrate_ode(%1%)", ts, "times",
89  static_cast<double*>(0));
90  stan::math::check_less("integrate_ode(%1%)",t0, ts[0], "initial time",
91  static_cast<double*>(0));
92 
93  const double absolute_tolerance = 1e-6;
94  const double relative_tolerance = 1e-6;
95  const double step_size = 0.1;
96 
97  // creates basic or coupled system by template specializations
99  coupled_system(f, y0, theta, x, x_int, msgs);
100 
101  // first time in the vector must be time of initial state
102  std::vector<double> ts_vec(ts.size() + 1);
103  ts_vec[0] = t0;
104  for (size_t n = 0; n < ts.size(); n++)
105  ts_vec[n+1] = ts[n];
106 
107  std::vector<std::vector<double> > y_coupled(ts_vec.size());
108  coupled_ode_observer observer(y_coupled);
109 
110  // the coupled system creates the coupled initial state
111  std::vector<double> initial_coupled_state
112  = coupled_system.initial_state();
113 
114  integrate_times(make_dense_output(absolute_tolerance,
115  relative_tolerance,
116  runge_kutta_dopri5<std::vector<double>,
117  double,
118  std::vector<double>,
119  double>() ),
120  coupled_system,
121  initial_coupled_state,
122  boost::begin(ts_vec), boost::end(ts_vec),
123  step_size,
124  observer);
125 
126  // remove the first state corresponding to the initial value
127  y_coupled.erase(y_coupled.begin());
128 
129  // the coupled system also encapsulates the decoupling operation
130  return coupled_system.decouple_states(y_coupled);
131  }
132 
133  }
134 
135 }
136 
137 #endif
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.
Definition: constants.hpp:86
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)
Definition: check_less.hpp:52
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.

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