Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
util.hpp
Go to the documentation of this file.
1 #ifndef STAN__MODEL__UTIL_HPP
2 #define STAN__MODEL__UTIL_HPP
3 
4 #include <cmath>
5 #include <iomanip>
6 #include <iostream>
7 #include <vector>
8 
11 #include <stan/agrad/rev/var.hpp>
12 #include <stan/agrad/rev/var.hpp>
13 #include <stan/agrad/autodiff.hpp>
14 
15 namespace stan {
16 
17  namespace model {
18 
40  template <bool jacobian_adjust_transform, class M>
41  double log_prob_propto(const M& model,
42  std::vector<double>& params_r,
43  std::vector<int>& params_i,
44  std::ostream* msgs = 0) {
45  using stan::agrad::var;
46  using std::vector;
47  vector<var> ad_params_r;
48  for (size_t i = 0; i < model.num_params_r(); ++i)
49  ad_params_r.push_back(params_r[i]);
50  try {
51  double lp
52  = model
53  .template log_prob<true,
54  jacobian_adjust_transform>(ad_params_r,params_i,
55  msgs)
56  .val();
58  return lp;
59  } catch (std::exception &ex) {
61  throw;
62  }
63  }
64 
82  template <bool propto, bool jacobian_adjust_transform, class M>
83  double log_prob_grad(const M& model,
84  std::vector<double>& params_r,
85  std::vector<int>& params_i,
86  std::vector<double>& gradient,
87  std::ostream* msgs = 0) {
88  using std::vector;
89  using stan::agrad::var;
90  double lp;
91  try {
92  vector<var> ad_params_r(params_r.size());
93  for (size_t i = 0; i < model.num_params_r(); ++i) {
94  stan::agrad::var var_i(params_r[i]);
95  ad_params_r[i] = var_i;
96  }
97  var adLogProb
98  = model
99  .template log_prob<propto,
100  jacobian_adjust_transform>(ad_params_r,
101  params_i,msgs);
102  lp = adLogProb.val();
103  adLogProb.grad(ad_params_r,gradient);
104  } catch (const std::exception &ex) {
106  throw;
107  }
109  return lp;
110  }
111 
132  template <bool jacobian_adjust_transform, class M>
133  double log_prob_propto(const M& model,
134  Eigen::VectorXd& params_r,
135  std::ostream* msgs = 0) {
136  using stan::agrad::var;
137  using std::vector;
138  double lp;
139  try {
140  vector<var> ad_params_r;
141  for (size_t i = 0; i < model.num_params_r(); ++i)
142  ad_params_r.push_back(params_r(i));
143  lp
144  = model
145  .template log_prob<true,
146  jacobian_adjust_transform>(ad_params_r, msgs)
147  .val();
148  } catch (std::exception &ex) {
150  throw;
151  }
153  return lp;
154  }
155 
172  template <bool propto, bool jacobian_adjust_transform, class M>
173  double log_prob_grad(const M& model,
174  Eigen::VectorXd& params_r,
175  Eigen::VectorXd& gradient,
176  std::ostream* msgs = 0) {
177  using std::vector;
178  using stan::agrad::var;
179  Eigen::Matrix<var,Eigen::Dynamic,1> ad_params_r(params_r.size());
180  for (size_t i = 0; i < model.num_params_r(); ++i) {
181  stan::agrad::var var_i(params_r[i]);
182  ad_params_r[i] = var_i;
183  }
184  try {
185  var adLogProb
186  = model
187  .template log_prob<propto,
188  jacobian_adjust_transform>(ad_params_r, msgs);
189  double val = adLogProb.val();
190  stan::agrad::grad(adLogProb, ad_params_r, gradient);
191  return val;
192  } catch (std::exception &ex) {
194  throw;
195  }
196  }
197 
216  template <bool propto, bool jacobian_adjust_transform, class M>
217  void finite_diff_grad(const M& model,
218  std::vector<double>& params_r,
219  std::vector<int>& params_i,
220  std::vector<double>& grad,
221  double epsilon = 1e-6,
222  std::ostream* msgs = 0) {
223  std::vector<double> perturbed(params_r);
224  grad.resize(params_r.size());
225  for (size_t k = 0; k < params_r.size(); k++) {
226  perturbed[k] += epsilon;
227  double logp_plus
228  = model
229  .template log_prob<propto,
230  jacobian_adjust_transform>(perturbed, params_i,
231  msgs);
232  perturbed[k] = params_r[k] - epsilon;
233  double logp_minus
234  = model
235  .template log_prob<propto,
236  jacobian_adjust_transform>(perturbed, params_i,
237  msgs);
238  double gradest = (logp_plus - logp_minus) / (2*epsilon);
239  grad[k] = gradest;
240  perturbed[k] = params_r[k];
241  }
242  }
243 
244 
267  template <bool propto, bool jacobian_adjust_transform, class M>
268  int test_gradients(const M& model,
269  std::vector<double>& params_r,
270  std::vector<int>& params_i,
271  double epsilon = 1e-6,
272  double error = 1e-6,
273  std::ostream& o = std::cout,
274  std::ostream* msgs = 0) {
275  std::vector<double> grad;
276  double lp
278  jacobian_adjust_transform>(model,
279  params_r,
280  params_i,
281  grad,msgs);
282 
283  std::vector<double> grad_fd;
284  finite_diff_grad<false,
285  true,
286  M>(model,
287  params_r, params_i,
288  grad_fd, epsilon,
289  msgs);
290 
291  int num_failed = 0;
292 
293  o << std::endl
294  << " Log probability=" << lp
295  << std::endl;
296 
297  o << std::endl
298  << std::setw(10) << "param idx"
299  << std::setw(16) << "value"
300  << std::setw(16) << "model"
301  << std::setw(16) << "finite diff"
302  << std::setw(16) << "error"
303  << std::endl;
304  for (size_t k = 0; k < params_r.size(); k++) {
305  o << std::setw(10) << k
306  << std::setw(16) << params_r[k]
307  << std::setw(16) << grad[k]
308  << std::setw(16) << grad_fd[k]
309  << std::setw(16) << (grad[k] - grad_fd[k])
310  << std::endl;
311  if (std::fabs(grad[k] - grad_fd[k]) > error)
312  num_failed++;
313  }
314  return num_failed;
315  }
316 
317 
340  template <bool propto, bool jacobian_adjust_transform, class M>
341  double grad_hess_log_prob(const M& model,
342  std::vector<double>& params_r,
343  std::vector<int>& params_i,
344  std::vector<double>& gradient,
345  std::vector<double>& hessian,
346  std::ostream* msgs = 0) {
347  const double epsilon = 1e-3;
348  const int order = 4;
349  const double perturbations[order]
350  = {-2*epsilon, -1*epsilon, epsilon, 2*epsilon};
351  const double coefficients[order]
352  = { 1.0 / 12.0,
353  -2.0 / 3.0,
354  2.0 / 3.0,
355  -1.0 / 12.0 };
356 
357  double result
359  jacobian_adjust_transform>(model,
360  params_r,
361  params_i,
362  gradient,
363  msgs);
364  hessian.assign(params_r.size() * params_r.size(), 0);
365  std::vector<double> temp_grad(params_r.size());
366  std::vector<double> perturbed_params(params_r.begin(), params_r.end());
367  for (size_t d = 0; d < params_r.size(); d++) {
368  double* row = &hessian[d*params_r.size()];
369  for (int i = 0; i < order; i++) {
370  perturbed_params[d] = params_r[d] + perturbations[i];
372  jacobian_adjust_transform>(model,
373  perturbed_params,
374  params_i,
375  temp_grad);
376  for (size_t dd = 0; dd < params_r.size(); dd++) {
377  row[dd] += 0.5 * coefficients[i] * temp_grad[dd] / epsilon;
378  hessian[d + dd*params_r.size()]
379  += 0.5 * coefficients[i] * temp_grad[dd] / epsilon;
380  }
381  }
382  perturbed_params[d] = params_r[d];
383  }
384  return result;
385  }
386 
387  // Interface for automatic differentiation of models
388 
389  template <class M>
391 
392  const M& model;
393  std::ostream* o;
394 
395  model_functional(const M& m, std::ostream* out): model(m), o(out) {};
396 
397  template <typename T>
398  T operator()(Eigen::Matrix<T, Eigen::Dynamic, 1>& x) const {
399  return model.template log_prob<true, true, T>(x, o);
400  }
401 
402  };
403 
404  template <class M>
405  void gradient(const M& model,
406  const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
407  double& f,
408  Eigen::Matrix<double, Eigen::Dynamic, 1>& grad_f,
409  std::ostream* msgs = 0) {
410 
411  stan::agrad::gradient(model_functional<M>(model, msgs), x, f, grad_f);
412 
413  }
414 
415  template <class M>
416  void hessian(const M& model,
417  const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
418  double& f,
419  Eigen::Matrix<double, Eigen::Dynamic, 1>& grad_f,
420  Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& hess_f,
421  std::ostream* msgs = 0) {
422 
423  stan::agrad::hessian(model_functional<M>(model, msgs), x, f, grad_f, hess_f);
424 
425  }
426 
427  template <class M>
428  void gradient_dot_vector(const M& model,
429  const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
430  const Eigen::Matrix<double, Eigen::Dynamic, 1>& v,
431  double& f,
432  double& grad_f_dot_v,
433  std::ostream* msgs = 0) {
434 
435  stan::agrad::gradient_dot_vector(model_functional<M>(model, msgs), x, v, f, grad_f_dot_v);
436 
437  }
438 
439  template <class M>
440  void hessian_times_vector(const M& model,
441  const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
442  const Eigen::Matrix<double, Eigen::Dynamic, 1>& v,
443  double& f,
444  Eigen::Matrix<double, Eigen::Dynamic, 1>& hess_f_dot_v,
445  std::ostream* msgs = 0) {
446 
447  stan::agrad::hessian_times_vector(model_functional<M>(model, msgs), x, v, f, hess_f_dot_v);
448 
449  }
450 
451  template <class M>
452  void grad_tr_mat_times_hessian(const M& model,
453  const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
454  const Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& X,
455  Eigen::Matrix<double, Eigen::Dynamic, 1>& grad_tr_X_hess_f,
456  std::ostream* msgs = 0) {
457 
458  stan::agrad::grad_tr_mat_times_hessian(model_functional<M>(model, msgs), x, X, grad_tr_X_hess_f);
459 
460  }
461 
462 
463  }
464 }
465 
466 
467 
468 #endif
void gradient(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &fx, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_fx)
Calculate the value and the gradient of the specified function at the specified argument.
Definition: autodiff.hpp:93
fvar< T > fabs(const fvar< T > &x)
Definition: fabs.hpp:18
static void recover_memory()
Recover memory used for all variables for reuse.
Definition: var_stack.hpp:52
double log_prob_grad(const M &model, std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &gradient, std::ostream *msgs=0)
Compute the gradient using reverse-mode automatic differentiation, writing the result into the specif...
Definition: util.hpp:83
Eigen::Matrix< T, 1, Eigen::Dynamic > row(const Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic > &m, size_t i)
Return the specified row of the specified matrix, using start-at-1 indexing.
Definition: row.hpp:26
void hessian(const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_f, Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &hess_f, std::ostream *msgs=0)
Definition: util.hpp:416
int test_gradients(const M &model, std::vector< double > &params_r, std::vector< int > &params_i, double epsilon=1e-6, double error=1e-6, std::ostream &o=std::cout, std::ostream *msgs=0)
Test the log_prob_grad() function's ability to produce accurate gradients using finite differences...
Definition: util.hpp:268
double val() const
Return the value of this variable.
Definition: var.hpp:209
static void grad(chainable *vi)
Compute the gradient for all variables starting from the specified root variable implementation.
Definition: chainable.hpp:110
void gradient(const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_f, std::ostream *msgs=0)
Definition: util.hpp:405
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)
Definition: autodiff.hpp:274
void hessian_times_vector(const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &v, double &f, Eigen::Matrix< double, Eigen::Dynamic, 1 > &hess_f_dot_v, std::ostream *msgs=0)
Definition: util.hpp:440
void finite_diff_grad(const M &model, std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &grad, double epsilon=1e-6, std::ostream *msgs=0)
Compute the gradient using finite differences for the specified parameters, writing the result into t...
Definition: util.hpp:217
model_functional(const M &m, std::ostream *out)
Definition: util.hpp:395
void grad_tr_mat_times_hessian(const F &f, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &M, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_tr_MH)
Definition: autodiff.hpp:321
void gradient_dot_vector(const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &v, double &f, double &grad_f_dot_v, std::ostream *msgs=0)
Definition: util.hpp:428
T operator()(Eigen::Matrix< T, Eigen::Dynamic, 1 > &x) const
Definition: util.hpp:398
void grad_tr_mat_times_hessian(const M &model, const Eigen::Matrix< double, Eigen::Dynamic, 1 > &x, const Eigen::Matrix< double, Eigen::Dynamic, Eigen::Dynamic > &X, Eigen::Matrix< double, Eigen::Dynamic, 1 > &grad_tr_X_hess_f, std::ostream *msgs=0)
Definition: util.hpp:452
double grad_hess_log_prob(const M &model, std::vector< double > &params_r, std::vector< int > &params_i, std::vector< double > &gradient, std::vector< double > &hessian, std::ostream *msgs=0)
Evaluate the log-probability, its gradient, and its Hessian at params_r.
Definition: util.hpp:341
Independent (input) and dependent (output) variables for gradients.
Definition: var.hpp:27
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:86
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)
Definition: autodiff.hpp:252
double log_prob_propto(const M &model, std::vector< double > &params_r, std::vector< int > &params_i, std::ostream *msgs=0)
Helper function to calculate log probability for double scalars up to a proportion.
Definition: util.hpp:41
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)
Definition: autodiff.hpp:194

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