-
-
Notifications
You must be signed in to change notification settings - Fork 182
/
hessian_times_vector.hpp
52 lines (46 loc) · 1.54 KB
/
hessian_times_vector.hpp
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
#ifndef STAN_MATH_MIX_FUNCTOR_HESSIAN_TIMES_VECTOR_HPP
#define STAN_MATH_MIX_FUNCTOR_HESSIAN_TIMES_VECTOR_HPP
#include <stan/math/fwd/core.hpp>
#include <stan/math/prim/fun/Eigen.hpp>
#include <stan/math/rev/core.hpp>
#include <stdexcept>
#include <vector>
namespace stan {
namespace math {
template <typename F>
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) {
using Eigen::Matrix;
// Run nested autodiff in this scope
nested_rev_autodiff nested;
Matrix<var, Eigen::Dynamic, 1> x_var(x.size());
for (int i = 0; i < x_var.size(); ++i) {
x_var(i) = x(i);
}
var fx_var;
var grad_fx_var_dot_v;
gradient_dot_vector(f, x_var, v, fx_var, grad_fx_var_dot_v);
fx = fx_var.val();
grad(grad_fx_var_dot_v.vi_);
Hv.resize(x.size());
for (int i = 0; i < x.size(); ++i) {
Hv(i) = x_var(i).adj();
}
}
template <typename T, typename F>
void hessian_times_vector(const F& f,
const Eigen::Matrix<T, Eigen::Dynamic, 1>& x,
const Eigen::Matrix<T, Eigen::Dynamic, 1>& v, T& fx,
Eigen::Matrix<T, Eigen::Dynamic, 1>& Hv) {
using Eigen::Matrix;
Matrix<T, Eigen::Dynamic, 1> grad;
Matrix<T, Eigen::Dynamic, Eigen::Dynamic> H;
hessian(f, x, fx, grad, H);
Hv = H * v;
}
} // namespace math
} // namespace stan
#endif