Stan  1.0
probability, sampling & optimization
newton.hpp
Go to the documentation of this file.
1 #ifndef __STAN__OPTIMIZATION__NEWTON_HPP__
2 #define __STAN__OPTIMIZATION__NEWTON_HPP__
3 
4 #include <Eigen/Dense>
5 #include <Eigen/Cholesky>
6 #include <Eigen/Eigenvalues>
8 
9 namespace stan {
10 
11  namespace optimization {
12 
13  typedef Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic> matrix_d;
14  typedef Eigen::Matrix<double, Eigen::Dynamic, 1> vector_d;
15 
16  namespace {
17  // Negates any positive eigenvalues in H so that H is negative
18  // definite, and then solves Hu = g and stores the result into
19  // g. Avoids problems due to non-log-concave distributions.
20  void make_negative_definite_and_solve(matrix_d& H, vector_d& g) {
21  Eigen::SelfAdjointEigenSolver<matrix_d> solver(H);
22  matrix_d eigenvectors = solver.eigenvectors();
23  vector_d eigenvalues = solver.eigenvalues();
24  vector_d eigenprojections = eigenvectors.transpose() * g;
25  for (int i = 0; i < g.size(); i++) {
26  eigenprojections[i] = -eigenprojections[i] / fabs(eigenvalues[i]);
27  }
28  g = eigenvectors * eigenprojections;
29  }
30  }
31 
33  std::vector<double>& params_r,
34  std::vector<int>& params_i,
35  std::ostream* output_stream = 0) {
36  std::vector<double> gradient;
37  std::vector<double> hessian;
38 
39  double f0 = model.grad_hess_log_prob(params_r, params_i,
40  gradient, hessian);
41  matrix_d H(params_r.size(), params_r.size());
42  for (size_t i = 0; i < hessian.size(); i++) {
43  H(i) = hessian[i];
44  }
45  vector_d g(params_r.size());
46  for (size_t i = 0; i < gradient.size(); i++)
47  g(i) = gradient[i];
48  make_negative_definite_and_solve(H, g);
49 // H.ldlt().solveInPlace(g);
50 
51  std::vector<double> new_params_r(params_r.size());
52  double step_size = 2;
53  double f1 = -1e100;
54  while (!(f1 >= f0)) {
55  step_size *= 0.5;
56  for (size_t i = 0; i < params_r.size(); i++)
57  new_params_r[i] = params_r[i] - step_size * g[i];
58  try {
59  f1 = model.grad_log_prob(new_params_r, params_i, gradient);
60  } catch (std::exception& e) {
61  f1 = -1e100;
62  }
63  }
64  for (size_t i = 0; i < params_r.size(); i++)
65  params_r[i] = new_params_r[i];
66 
67  return f1;
68  }
69 
70  }
71 
72 }
73 
74 #endif
The prob_grad class represents densities with fixed numbers of discrete and scalar parameters and the...
Definition: prob_grad.hpp:24
virtual double grad_log_prob(std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &gradient, std::ostream *output_stream=0)=0
virtual double grad_hess_log_prob(std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &gradient, std::vector< double > &hessian, std::ostream *output_stream=0)
Evaluate the log-probability, its gradient, and its Hessian at params_r.
Definition: prob_grad.hpp:112
var fabs(const var &a)
Return the absolute value of the variable (cmath).
Definition: agrad.hpp:2023
double e()
Return the base of the natural logarithm.
double newton_step(stan::model::prob_grad &model, std::vector< double > &params_r, std::vector< int > &params_i, std::ostream *output_stream=0)
Definition: newton.hpp:32
Eigen::Matrix< double, Eigen::Dynamic, 1 > vector_d
Definition: newton.hpp:14
Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > matrix_d
Definition: newton.hpp:13
Probability, optimization and sampling library.
Definition: agrad.cpp:6

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