Stan Math Library  2.10.0
reverse mode automatic differentiation
hessian.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_MIX_MAT_FUNCTOR_HESSIAN_HPP
2 #define STAN_MATH_MIX_MAT_FUNCTOR_HESSIAN_HPP
3 
4 #include <stan/math/fwd/core.hpp>
6 #include <stan/math/rev/core.hpp>
7 #include <stdexcept>
8 #include <vector>
9 
10 namespace stan {
11 
12  namespace math {
13 
43  template <typename F>
44  void
45  hessian(const F& f,
46  const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
47  double& fx,
48  Eigen::Matrix<double, Eigen::Dynamic, 1>& grad,
49  Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& H) {
50  H.resize(x.size(), x.size());
51  grad.resize(x.size());
52  try {
53  for (int i = 0; i < x.size(); ++i) {
54  start_nested();
55  Eigen::Matrix<fvar<var>, Eigen::Dynamic, 1> x_fvar(x.size());
56  for (int j = 0; j < x.size(); ++j)
57  x_fvar(j) = fvar<var>(x(j), i == j);
58  fvar<var> fx_fvar = f(x_fvar);
59  grad(i) = fx_fvar.d_.val();
60  if (i == 0) fx = fx_fvar.val_.val();
61  stan::math::grad(fx_fvar.d_.vi_);
62  for (int j = 0; j < x.size(); ++j)
63  H(i, j) = x_fvar(j).val_.adj();
65  }
66  } catch (const std::exception& e) {
68  throw;
69  }
70  }
71  // time O(N^3); space O(N^2)
72  template <typename T, typename F>
73  void
74  hessian(const F& f,
75  const Eigen::Matrix<T, Eigen::Dynamic, 1>& x,
76  T& fx,
77  Eigen::Matrix<T, Eigen::Dynamic, 1>& grad,
78  Eigen::Matrix<T, Eigen::Dynamic, Eigen::Dynamic>& H) {
79  H.resize(x.size(), x.size());
80  grad.resize(x.size());
81  Eigen::Matrix<fvar<fvar<T> >, Eigen::Dynamic, 1> x_fvar(x.size());
82  for (int i = 0; i < x.size(); ++i) {
83  for (int j = i; j < x.size(); ++j) {
84  for (int k = 0; k < x.size(); ++k)
85  x_fvar(k) = fvar<fvar<T> >(fvar<T>(x(k), j == k),
86  fvar<T>(i == k, 0));
87  fvar<fvar<T> > fx_fvar = f(x_fvar);
88  if (j == 0)
89  fx = fx_fvar.val_.val_;
90  if (i == j)
91  grad(i) = fx_fvar.d_.val_;
92  H(i, j) = fx_fvar.d_.d_;
93  H(j, i) = H(i, j);
94  }
95  }
96  }
97 
98  } // namespace math
99 } // namespace stan
100 #endif
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)
Calculate the value, the gradient, and the Hessian, of the specified function at the specified argume...
Definition: hessian.hpp:45
static void grad(vari *vi)
Compute the gradient for all variables starting from the specified root variable implementation.
Definition: grad.hpp:30
vari * vi_
Pointer to the implementation of this variable.
Definition: var.hpp:43
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:95
static void recover_memory_nested()
Recover only the memory used for the top nested call.
double val() const
Return the value of this variable.
Definition: var.hpp:233
static void start_nested()
Record the current position so that recover_memory_nested() can find it.

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