Stan Math Library  2.10.0
reverse mode automatic differentiation
hessian_times_vector.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_MIX_MAT_FUNCTOR_HESSIAN_TIMES_VECTOR_HPP
2 #define STAN_MATH_MIX_MAT_FUNCTOR_HESSIAN_TIMES_VECTOR_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 
14  template <typename F>
15  void
17  const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
18  const Eigen::Matrix<double, Eigen::Dynamic, 1>& v,
19  double& fx,
20  Eigen::Matrix<double, Eigen::Dynamic, 1>& Hv) {
21  using stan::math::fvar;
22  using stan::math::var;
23  using Eigen::Matrix;
24  start_nested();
25  try {
26  Matrix<var, Eigen::Dynamic, 1> x_var(x.size());
27  for (int i = 0; i < x_var.size(); ++i)
28  x_var(i) = x(i);
29  var fx_var;
30  var grad_fx_var_dot_v;
31  gradient_dot_vector(f, x_var, v, fx_var, grad_fx_var_dot_v);
32  fx = fx_var.val();
33  stan::math::grad(grad_fx_var_dot_v.vi_);
34  Hv.resize(x.size());
35  for (int i = 0; i < x.size(); ++i)
36  Hv(i) = x_var(i).adj();
37  } catch (const std::exception& e) {
39  throw;
40  }
42  }
43  template <typename T, typename F>
44  void
46  const Eigen::Matrix<T, Eigen::Dynamic, 1>& x,
47  const Eigen::Matrix<T, Eigen::Dynamic, 1>& v,
48  T& fx,
49  Eigen::Matrix<T, Eigen::Dynamic, 1>& Hv) {
50  using Eigen::Matrix;
51  Matrix<T, Eigen::Dynamic, 1> grad;
52  Matrix<T, Eigen::Dynamic, Eigen::Dynamic> H;
53  hessian(f, x, fx, grad, H);
54  Hv = H * v;
55  }
56 
57  } // namespace math
58 } // namespace stan
59 #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
Independent (input) and dependent (output) variables for gradients.
Definition: var.hpp:31
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)
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
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)
static void start_nested()
Record the current position so that recover_memory_nested() can find it.

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