Stan Math Library  2.20.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 namespace math {
12 
13 template <typename F>
14 void hessian_times_vector(const F& f,
15  const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
16  const Eigen::Matrix<double, Eigen::Dynamic, 1>& v,
17  double& fx,
18  Eigen::Matrix<double, Eigen::Dynamic, 1>& Hv) {
19  using Eigen::Matrix;
20  start_nested();
21  try {
22  Matrix<var, Eigen::Dynamic, 1> x_var(x.size());
23  for (int i = 0; i < x_var.size(); ++i)
24  x_var(i) = x(i);
25  var fx_var;
26  var grad_fx_var_dot_v;
27  gradient_dot_vector(f, x_var, v, fx_var, grad_fx_var_dot_v);
28  fx = fx_var.val();
29  grad(grad_fx_var_dot_v.vi_);
30  Hv.resize(x.size());
31  for (int i = 0; i < x.size(); ++i)
32  Hv(i) = x_var(i).adj();
33  } catch (const std::exception& e) {
35  throw;
36  }
38 }
39 template <typename T, typename F>
40 void hessian_times_vector(const F& f,
41  const Eigen::Matrix<T, Eigen::Dynamic, 1>& x,
42  const Eigen::Matrix<T, Eigen::Dynamic, 1>& v, T& fx,
43  Eigen::Matrix<T, Eigen::Dynamic, 1>& Hv) {
44  using Eigen::Matrix;
45  Matrix<T, Eigen::Dynamic, 1> grad;
46  Matrix<T, Eigen::Dynamic, Eigen::Dynamic> H;
47  hessian(f, x, fx, grad, H);
48  Hv = H * v;
49 }
50 
51 } // namespace math
52 } // namespace stan
53 #endif
void hessian(const F &f, const Eigen::Matrix< T, Eigen::Dynamic, 1 > &x, T &fx, Eigen::Matrix< T, Eigen::Dynamic, 1 > &grad, Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &H)
Calculate the value, the gradient, and the Hessian, of the specified function at the specified argume...
Definition: hessian.hpp:41
Independent (input) and dependent (output) variables for gradients.
Definition: var.hpp:33
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:45
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:87
static void recover_memory_nested()
Recover only the memory used for the top nested call.
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.
double val() const
Return the value of this variable.
Definition: var.hpp:294

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