Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
coupled_ode_system.hpp
Go to the documentation of this file.
1 #ifndef STAN__MATH__ODE__COUPLED_ODE_SYSTEM_HPP
2 #define STAN__MATH__ODE__COUPLED_ODE_SYSTEM_HPP
3 
4 #include <ostream>
5 #include <vector>
6 
8 
9 namespace stan {
10 
11  namespace math {
12 
26  template <typename F, typename T1, typename T2>
28  };
29 
30 
41  template <typename F>
42  struct coupled_ode_system<F, double, double> {
43 
44  const F& f_;
45  const std::vector<double>& y0_dbl_;
46  const std::vector<double>& theta_dbl_;
47  const std::vector<double>& x_;
48  const std::vector<int>& x_int_;
49  const int N_;
50  const int M_;
51  const int size_;
52  std::ostream* msgs_;
53 
66  coupled_ode_system(const F& f,
67  const std::vector<double>& y0,
68  const std::vector<double>& theta,
69  const std::vector<double>& x,
70  const std::vector<int>& x_int,
71  std::ostream* msgs)
72  : f_(f),
73  y0_dbl_(y0),
74  theta_dbl_(theta),
75  x_(x),
76  x_int_(x_int),
77  N_(y0.size()),
78  M_(theta.size()),
79  size_(N_),
80  msgs_(msgs) {
81  }
82 
98  void operator()(const std::vector<double>& y,
99  std::vector<double>& dy_dt,
100  double t) {
101  dy_dt = f_(t,y,theta_dbl_,x_,x_int_,msgs_);
102  stan::math::check_matching_sizes("coupled_ode_system(%1%)",
103  y,"y",dy_dt,"dy_dt",
104  static_cast<double*>(0));
105  }
106 
112  int size() const {
113  return size_;
114  }
115 
128  std::vector<double> initial_state() {
129  std::vector<double> state(size_, 0.0);
130  for (int n = 0; n < N_; n++)
131  state[n] = y0_dbl_[n];
132  return state;
133  }
134 
145  std::vector<std::vector<double> >
146  decouple_states(const std::vector<std::vector<double> >& y) {
147  return y;
148  }
149 
150  };
151 
152  }
153 
154 }
155 
156 #endif
int N_
void operator()(const std::vector< double > &y, std::vector< double > &dy_dt, double t)
Calculates the derivative of the coupled ode system with respect to the specified state at the specif...
std::vector< double > initial_state()
Returns the initial state of the coupled system, which is identical to the base ODE original state in...
int size() const
Returns the size of the coupled system.
coupled_ode_system(const F &f, const std::vector< double > &y0, const std::vector< double > &theta, const std::vector< double > &x, const std::vector< int > &x_int, std::ostream *msgs)
Construct the coupled ODE system from the base system function, initial state, parameters, data and a stream for messages.
std::vector< std::vector< double > > decouple_states(const std::vector< std::vector< double > > &y)
Returns the base portion of the coupled state.
int size(const std::vector< T > &x)
Definition: size.hpp:11
size_t size_
Definition: dot_self.hpp:18
Base template class for a coupled ordinary differential equation system, which adds sensitivities to ...
bool check_matching_sizes(const char *function, const T_y1 &y1, const char *name1, const T_y2 &y2, const char *name2, T_result *result)
int M_

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