Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
autodiff.hpp
Go to the documentation of this file.
1 #ifndef STAN__AGRAD__AUTO_DIFF_HPP
2 #define STAN__AGRAD__AUTO_DIFF_HPP
3 
5 #include <stan/agrad/rev.hpp>
6 #include <stan/agrad/fwd.hpp>
7 
8 namespace stan {
9 
10  namespace agrad {
11 
23  template <typename T, typename F>
24  void
25  derivative(const F& f,
26  const T& x,
27  T& fx,
28  T& dfx_dx) {
29  fvar<T> x_fvar = fvar<T>(x,1.0);
30  fvar<T> fx_fvar = f(x_fvar);
31  fx = fx_fvar.val_;
32  dfx_dx = fx_fvar.d_;
33  }
34 
47  template <typename T, typename F>
48  void
49  partial_derivative(const F& f,
50  const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
51  int n,
52  T& fx,
53  T& dfx_dxn) {
54  Eigen::Matrix<fvar<T>,Eigen::Dynamic,1> x_fvar(x.size());
55  for (int i = 0; i < x.size(); ++i)
56  x_fvar(i) = fvar<T>(x(i),i==n);
57  fvar<T> fx_fvar = f(x_fvar);
58  fx = fx_fvar.val_;
59  dfx_dxn = fx_fvar.d_;
60  }
61 
91  template <typename F>
92  void
93  gradient(const F& f,
94  const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
95  double& fx,
96  Eigen::Matrix<double,Eigen::Dynamic,1>& grad_fx) {
97  using stan::agrad::var;
98  start_nested();
99  try {
100  Eigen::Matrix<var,Eigen::Dynamic,1> x_var(x.size());
101  for (int i = 0; i < x.size(); ++i)
102  x_var(i) = x(i);
103  var fx_var = f(x_var);
104  fx = fx_var.val();
105  grad_fx.resize(x.size());
106  stan::agrad::grad(fx_var.vi_);
107  for (int i = 0; i < x.size(); ++i)
108  grad_fx(i) = x_var(i).adj();
109  } catch (const std::exception& /*e*/) {
111  throw;
112  }
114  }
115  template <typename T, typename F>
116  void
117  gradient(const F& f,
118  const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
119  T& fx,
120  Eigen::Matrix<T,Eigen::Dynamic,1>& grad_fx) {
121  Eigen::Matrix<fvar<T>,Eigen::Dynamic,1> x_fvar(x.size());
122  grad_fx.resize(x.size());
123  for (int i = 0; i < x.size(); ++i) {
124  for (int k = 0; k < x.size(); ++k)
125  x_fvar(k) = fvar<T>(x(k),k==i);
126  fvar<T> fx_fvar = f(x_fvar);
127  if (i == 0) fx = fx_fvar.val_;
128  grad_fx(i) = fx_fvar.d_;
129  }
130  }
131 
132 
133  template <typename F>
134  void
135  jacobian(const F& f,
136  const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
137  Eigen::Matrix<double,Eigen::Dynamic,1>& fx,
138  Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& J) {
139  using Eigen::Matrix; using Eigen::Dynamic;
140  using stan::agrad::var;
141  start_nested();
142  try {
143  Matrix<var,Dynamic,1> x_var(x.size());
144  for (int k = 0; k < x.size(); ++k)
145  x_var(k) = x(k);
146  Matrix<var,Dynamic,1> fx_var = f(x_var);
147  fx.resize(fx_var.size());
148  for (int i = 0; i < fx_var.size(); ++i)
149  fx(i) = fx_var(i).val();
150  J.resize(x.size(), fx_var.size());
151  for (int i = 0; i < fx_var.size(); ++i) {
152  if (i > 0)
154  grad(fx_var(i).vi_);
155  for (int k = 0; k < x.size(); ++k)
156  J(k,i) = x_var(k).adj();
157  }
158  } catch (const std::exception& e) {
160  throw;
161  }
163  }
164  template <typename T, typename F>
165  void
166  jacobian(const F& f,
167  const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
168  Eigen::Matrix<T,Eigen::Dynamic,1>& fx,
169  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& J) {
170  using Eigen::Matrix; using Eigen::Dynamic;
171  using stan::agrad::fvar;
172  Matrix<fvar<T>,Dynamic,1> x_fvar(x.size());
173  for (int i = 0; i < x.size(); ++i) {
174  for (int k = 0; k < x.size(); ++k)
175  x_fvar(k) = fvar<T>(x(k), i == k);
176  Matrix<fvar<T>,Dynamic,1> fx_fvar
177  = f(x_fvar);
178  if (i == 0) {
179  J.resize(x.size(),fx_fvar.size());
180  fx.resize(fx_fvar.size());
181  for (int k = 0; k < fx_fvar.size(); ++k)
182  fx(k) = fx_fvar(k).val_;
183  }
184  for (int k = 0; k < fx_fvar.size(); ++k) {
185  J(i,k) = fx_fvar(k).d_;
186  }
187  }
188  }
189 
190 
191  // time O(N^2); space O(N^2)
192  template <typename F>
193  void
194  hessian(const F& f,
195  const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
196  double& fx,
197  Eigen::Matrix<double,Eigen::Dynamic,1>& grad,
198  Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& H) {
199  start_nested();
200  try {
201  H.resize(x.size(), x.size());
202  grad.resize(x.size());
203  for (int i = 0; i < x.size(); ++i) {
204  Eigen::Matrix<fvar<var>, Eigen::Dynamic, 1> x_fvar(x.size());
205  for (int j = 0; j < x.size(); ++j)
206  x_fvar(j) = fvar<var>(x(j),i==j);
207  fvar<var> fx_fvar = f(x_fvar);
208  grad(i) = fx_fvar.d_.val();
209  if (i == 0) fx = fx_fvar.val_.val();
210  stan::agrad::grad(fx_fvar.d_.vi_);
211  for (int j = 0; j < x.size(); ++j)
212  H(i,j) = x_fvar(j).val_.adj();
213  }
214  } catch (const std::exception& e) {
216  throw;
217  }
219  }
220  // time O(N^3); space O(N^2)
221  template <typename T, typename F>
222  void
223  hessian(const F& f,
224  const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
225  T& fx,
226  Eigen::Matrix<T,Eigen::Dynamic,1>& grad,
227  Eigen::Matrix<T,Eigen::Dynamic,Eigen::Dynamic>& H) {
228  H.resize(x.size(), x.size());
229  grad.resize(x.size());
230  Eigen::Matrix<fvar<fvar<T> >,Eigen::Dynamic,1> x_fvar(x.size());
231  for (int i = 0; i < x.size(); ++i) {
232  for (int j = i; j < x.size(); ++j) {
233  for (int k = 0; k < x.size(); ++k)
234  x_fvar(k) = fvar<fvar<T> >(fvar<T>(x(k),j==k),
235  fvar<T>(i==k,0));
236  fvar<fvar<T> > fx_fvar = f(x_fvar);
237  if (j == 0)
238  fx = fx_fvar.val_.val_;
239  if (i == j)
240  grad(i) = fx_fvar.d_.val_;
241  H(i,j) = fx_fvar.d_.d_;
242  H(j,i) = H(i,j);
243  }
244  }
245  }
246 
247 
248  // aka directional derivative (not length normalized)
249  // T2 must be assignable to T1
250  template <typename T1, typename T2, typename F>
251  void
253  const Eigen::Matrix<T1,Eigen::Dynamic,1>& x,
254  const Eigen::Matrix<T2,Eigen::Dynamic,1>& v,
255  T1& fx,
256  T1& grad_fx_dot_v) {
257  using stan::agrad::fvar;
258  using stan::agrad::var;
259  using Eigen::Matrix;
260  using Eigen::Dynamic;
261  Matrix<fvar<T1>,Dynamic,1> x_fvar(x.size());
262  for (int i = 0; i < x.size(); ++i)
263  x_fvar(i) = fvar<T1>(x(i),v(i));
264  fvar<T1> fx_fvar = f(x_fvar);
265  fx = fx_fvar.val_;
266  grad_fx_dot_v = fx_fvar.d_;
267  }
268 
269 
270 
271 
272  template <typename F>
273  void
275  const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
276  const Eigen::Matrix<double,Eigen::Dynamic,1>& v,
277  double& fx,
278  Eigen::Matrix<double,Eigen::Dynamic,1>& Hv) {
279  using stan::agrad::fvar;
280  using stan::agrad::var;
281  using Eigen::Matrix;
282  using Eigen::Dynamic;
283  start_nested();
284  try {
285  Matrix<var,Dynamic,1> x_var(x.size());
286  for (int i = 0; i < x_var.size(); ++i)
287  x_var(i) = x(i);
288  var fx_var;
289  var grad_fx_var_dot_v;
290  gradient_dot_vector(f,x_var,v,fx_var,grad_fx_var_dot_v);
291  fx = fx_var.val();
292  stan::agrad::grad(grad_fx_var_dot_v.vi_);
293  Hv.resize(x.size());
294  for (int i = 0; i < x.size(); ++i)
295  Hv(i) = x_var(i).adj();
296  } catch (const std::exception& e) {
298  throw;
299  }
301  }
302  template <typename T, typename F>
303  void
305  const Eigen::Matrix<T,Eigen::Dynamic,1>& x,
306  const Eigen::Matrix<T,Eigen::Dynamic,1>& v,
307  T& fx,
308  Eigen::Matrix<T,Eigen::Dynamic,1>& Hv) {
309  using Eigen::Matrix;
310  using Eigen::Dynamic;
311  Matrix<T,Dynamic,1> grad;
312  Matrix<T,Dynamic,Dynamic> H;
313  hessian(f,x,fx,grad,H);
314  Hv = H * v;
315  }
316 
317  // FIXME: add other results that are easy to extract
318  // // N * (fwd(2) + bk)
319  template <typename F>
320  void
322  const F& f,
323  const Eigen::Matrix<double,Eigen::Dynamic,1>& x,
324  const Eigen::Matrix<double,Eigen::Dynamic,Eigen::Dynamic>& M,
325  Eigen::Matrix<double,Eigen::Dynamic,1>& grad_tr_MH) {
326  using Eigen::Matrix;
327  using Eigen::Dynamic;
328  start_nested();
329  try {
330 
331  grad_tr_MH.resize(x.size());
332 
333  Matrix<var,Dynamic,1> x_var(x.size());
334  for (int i = 0; i < x.size(); ++i)
335  x_var(i) = x(i);
336 
337  Matrix<fvar<var>,Dynamic,1> x_fvar(x.size());
338 
339  var sum(0.0);
340  Matrix<double,Dynamic,1> M_n(x.size());
341  for (int n = 0; n < x.size(); ++n) {
342  for (int k = 0; k < x.size(); ++k)
343  M_n(k) = M(n,k);
344  for (int k = 0; k < x.size(); ++k)
345  x_fvar(k) = fvar<var>(x_var(k), k == n);
346  fvar<var> fx;
347  fvar<var> grad_fx_dot_v;
348  gradient_dot_vector<fvar<var>,double>(f,x_fvar,M_n,fx,grad_fx_dot_v);
349  sum += grad_fx_dot_v.d_;
350  }
351 
352  stan::agrad::grad(sum.vi_);
353  for (int i = 0; i < x.size(); ++i)
354  grad_tr_MH(i) = x_var(i).adj();
355  } catch (const std::exception& e) {
357  throw;
358  }
360  }
361 
362  }
363 }
364 #endif
void gradient(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_fx)
Calculate the value and the gradient of the specified function at the specified argument.
Definition: autodiff.hpp:93
static void set_zero_all_adjoints()
Reset all adjoint values in the stack to zero.
Definition: chainable.hpp:87
static void start_nested()
Record the current position so that recover_memory_nested() can find it.
Definition: var_stack.hpp:96
double val() const
Return the value of this variable.
Definition: var.hpp:209
static void grad(chainable *vi)
Compute the gradient for all variables starting from the specified root variable implementation.
Definition: chainable.hpp:110
void partial_derivative(const F &f, const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, int n, T &fx, T &dfx_dxn)
Return the partial derivative of the specified multiivariate function at the specified argument...
Definition: autodiff.hpp:49
void hessian_times_vector(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &v, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &Hv)
Definition: autodiff.hpp:274
fvar< T > sum(const Eigen::Matrix< fvar< T >, R, C > &m)
Definition: sum.hpp:14
vari * vi_
Pointer to the implementation of this variable.
Definition: var.hpp:40
void derivative(const F &f, const T &x, T &fx, T &dfx_dx)
Return the derivative of the specified univariate function at the specified argument.
Definition: autodiff.hpp:25
void grad_tr_mat_times_hessian(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &M, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_tr_MH)
Definition: autodiff.hpp:321
Independent (input) and dependent (output) variables for gradients.
Definition: var.hpp:27
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:86
void jacobian(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, Eigen::Matrix< double, Eigen::Dynamic, 1 > &fx, Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &J)
Definition: autodiff.hpp:135
void gradient_dot_vector(const F &f, const Eigen::Matrix< T1, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< T2, Eigen::Dynamic, 1 > &v, T1 &fx, T1 &grad_fx_dot_v)
Definition: autodiff.hpp:252
static void recover_memory_nested()
Recover only the memory used for the top nested call.
Definition: var_stack.hpp:72
void hessian(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad, Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &H)
Definition: autodiff.hpp:194

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