1 #ifndef STAN__MODEL__UTIL_HPP
2 #define STAN__MODEL__UTIL_HPP
40 template <
bool jacobian_adjust_transform,
class M>
42 std::vector<double>& params_r,
43 std::vector<int>& params_i,
44 std::ostream* msgs = 0) {
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]);
53 .template log_prob<
true,
54 jacobian_adjust_transform>(ad_params_r,params_i,
59 }
catch (std::exception &ex) {
82 template <
bool propto,
bool jacobian_adjust_transform,
class M>
84 std::vector<double>& params_r,
85 std::vector<int>& params_i,
87 std::ostream* msgs = 0) {
92 vector<var> ad_params_r(params_r.size());
93 for (
size_t i = 0; i < model.num_params_r(); ++i) {
95 ad_params_r[i] = var_i;
99 .template log_prob<propto,
100 jacobian_adjust_transform>(ad_params_r,
102 lp = adLogProb.
val();
103 adLogProb.grad(ad_params_r,gradient);
104 }
catch (
const std::exception &ex) {
132 template <
bool jacobian_adjust_transform,
class M>
134 Eigen::VectorXd& params_r,
135 std::ostream* msgs = 0) {
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));
145 .template log_prob<
true,
146 jacobian_adjust_transform>(ad_params_r, msgs)
148 }
catch (std::exception &ex) {
172 template <
bool propto,
bool jacobian_adjust_transform,
class M>
174 Eigen::VectorXd& params_r,
176 std::ostream* msgs = 0) {
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) {
182 ad_params_r[i] = var_i;
187 .template log_prob<propto,
188 jacobian_adjust_transform>(ad_params_r, msgs);
189 double val = adLogProb.
val();
192 }
catch (std::exception &ex) {
216 template <
bool propto,
bool jacobian_adjust_transform,
class M>
218 std::vector<double>& params_r,
219 std::vector<int>& params_i,
220 std::vector<double>&
grad,
221 double epsilon = 1
e-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;
229 .template log_prob<propto,
230 jacobian_adjust_transform>(perturbed, params_i,
232 perturbed[k] = params_r[k] - epsilon;
235 .template log_prob<propto,
236 jacobian_adjust_transform>(perturbed, params_i,
238 double gradest = (logp_plus - logp_minus) / (2*epsilon);
240 perturbed[k] = params_r[k];
267 template <
bool propto,
bool jacobian_adjust_transform,
class M>
269 std::vector<double>& params_r,
270 std::vector<int>& params_i,
271 double epsilon = 1
e-6,
273 std::ostream& o = std::cout,
274 std::ostream* msgs = 0) {
275 std::vector<double>
grad;
278 jacobian_adjust_transform>(model,
283 std::vector<double> grad_fd;
294 <<
" Log probability=" << lp
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"
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])
311 if (
std::fabs(grad[k] - grad_fd[k]) > error)
340 template <
bool propto,
bool jacobian_adjust_transform,
class M>
342 std::vector<double>& params_r,
343 std::vector<int>& params_i,
346 std::ostream* msgs = 0) {
347 const double epsilon = 1
e-3;
349 const double perturbations[order]
350 = {-2*epsilon, -1*epsilon, epsilon, 2*epsilon};
351 const double coefficients[order]
359 jacobian_adjust_transform>(model,
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,
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;
382 perturbed_params[d] = params_r[d];
397 template <
typename T>
399 return model.template log_prob<true, true, T>(x,
o);
406 const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
408 Eigen::Matrix<double, Eigen::Dynamic, 1>& grad_f,
409 std::ostream* msgs = 0) {
417 const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
419 Eigen::Matrix<double, Eigen::Dynamic, 1>& grad_f,
420 Eigen::Matrix<double, Eigen::Dynamic, Eigen::Dynamic>& hess_f,
421 std::ostream* msgs = 0) {
429 const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
430 const Eigen::Matrix<double, Eigen::Dynamic, 1>& v,
432 double& grad_f_dot_v,
433 std::ostream* msgs = 0) {
441 const Eigen::Matrix<double, Eigen::Dynamic, 1>& x,
442 const Eigen::Matrix<double, Eigen::Dynamic, 1>& v,
444 Eigen::Matrix<double, Eigen::Dynamic, 1>& hess_f_dot_v,
445 std::ostream* msgs = 0) {
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) {
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.
fvar< T > fabs(const fvar< T > &x)
static void recover_memory()
Recover memory used for all variables for reuse.
double log_prob_grad(const M &model, std::vector< double > ¶ms_r, std::vector< int > ¶ms_i, std::vector< double > &gradient, std::ostream *msgs=0)
Compute the gradient using reverse-mode automatic differentiation, writing the result into the specif...
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.
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)
int test_gradients(const M &model, std::vector< double > ¶ms_r, std::vector< int > ¶ms_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...
double val() const
Return the value of this variable.
static void grad(chainable *vi)
Compute the gradient for all variables starting from the specified root variable implementation.
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)
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)
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)
void finite_diff_grad(const M &model, std::vector< double > ¶ms_r, std::vector< int > ¶ms_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...
model_functional(const M &m, std::ostream *out)
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)
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)
T operator()(Eigen::Matrix< T, Eigen::Dynamic, 1 > &x) const
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)
double grad_hess_log_prob(const M &model, std::vector< double > ¶ms_r, std::vector< int > ¶ms_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.
Independent (input) and dependent (output) variables for gradients.
double e()
Return the base of the natural logarithm.
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)
double log_prob_propto(const M &model, std::vector< double > ¶ms_r, std::vector< int > ¶ms_i, std::ostream *msgs=0)
Helper function to calculate log probability for double scalars up to a proportion.
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)