Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
Classes | Functions
stan::model Namespace Reference

For compiling models. More...

Classes

class  prob_grad
 The prob_grad class represents the basic parameter holders for a model. More...
 
struct  model_functional
 

Functions

template<bool jacobian_adjust_transform, class M >
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. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
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 specified gradient, using the specified perturbation. More...
 
template<bool jacobian_adjust_transform, class M >
double log_prob_propto (const M &model, Eigen::VectorXd &params_r, std::ostream *msgs=0)
 Helper function to calculate log probability for double scalars up to a proportion. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
double log_prob_grad (const M &model, Eigen::VectorXd &params_r, Eigen::VectorXd &gradient, std::ostream *msgs=0)
 Compute the gradient using reverse-mode automatic differentiation, writing the result into the specified gradient, using the specified perturbation. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
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 the specified gradient, using the specified perturbation. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
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. More...
 
template<bool propto, bool jacobian_adjust_transform, class M >
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. More...
 
template<class M >
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)
 
template<class M >
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)
 
template<class M >
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)
 
template<class M >
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)
 
template<class M >
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)
 

Detailed Description

For compiling models.

Function Documentation

template<bool propto, bool jacobian_adjust_transform, class M >
void stan::model::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 the specified gradient, using the specified perturbation.

Template Parameters
proptoTrue if calculation is up to proportion (double-only terms dropped).
jacobian_adjust_transformTrue if the log absolute Jacobian determinant of inverse parameter transforms is added to the log probability.
MClass of model.
Parameters
modelModel.
params_rReal-valued parameters.
params_iInteger-valued parameters.
[out]gradVector into which gradient is written.
epsilon
[in,out]msgs

Definition at line 217 of file util.hpp.

template<bool propto, bool jacobian_adjust_transform, class M >
double stan::model::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.

This default version computes the Hessian numerically by finite-differencing the gradient, at a cost of O(params_r.size()^2).

Template Parameters
proptoTrue if calculation is up to proportion (double-only terms dropped).
jacobian_adjust_transformTrue if the log absolute Jacobian determinant of inverse parameter transforms is added to the log probability.
MClass of model.
Parameters
modelModel.
params_rReal-valued parameter vector.
params_iInteger-valued parameter vector.
gradientVector to write gradient to.
hessianVector to write gradient to. hessian[i*D + j] gives the element at the ith row and jth column of the Hessian (where D=params_r.size()).
msgsStream to which print statements in Stan programs are written, default is 0

Definition at line 341 of file util.hpp.

template<class M >
void stan::model::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 at line 452 of file util.hpp.

template<class M >
void stan::model::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 at line 405 of file util.hpp.

template<class M >
void stan::model::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 at line 428 of file util.hpp.

template<class M >
void stan::model::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 at line 416 of file util.hpp.

template<class M >
void stan::model::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 at line 440 of file util.hpp.

template<bool propto, bool jacobian_adjust_transform, class M >
double stan::model::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 specified gradient, using the specified perturbation.

Template Parameters
proptoTrue if calculation is up to proportion (double-only terms dropped).
jacobian_adjust_transformTrue if the log absolute Jacobian determinant of inverse parameter transforms is added to the log probability.
MClass of model.
Parameters
[in]modelModel.
[in]params_rReal-valued parameters.
[in]params_iInteger-valued parameters.
[out]gradientVector into which gradient is written.
[in,out]msgs

Definition at line 83 of file util.hpp.

template<bool propto, bool jacobian_adjust_transform, class M >
double stan::model::log_prob_grad ( const M &  model,
Eigen::VectorXd &  params_r,
Eigen::VectorXd &  gradient,
std::ostream *  msgs = 0 
)

Compute the gradient using reverse-mode automatic differentiation, writing the result into the specified gradient, using the specified perturbation.

Template Parameters
proptoTrue if calculation is up to proportion (double-only terms dropped).
jacobian_adjust_transformTrue if the log absolute Jacobian determinant of inverse parameter transforms is added to the log probability.
MClass of model.
Parameters
[in]modelModel.
[in]params_rReal-valued parameters.
[out]gradientVector into which gradient is written.
[in,out]msgs

Definition at line 173 of file util.hpp.

template<bool jacobian_adjust_transform, class M >
double stan::model::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.

This implementation wraps the double values in stan::agrad::var and calls the model's log_prob() function with propto=true and the specified parameter for applying the Jacobian adjustment for transformed parameters.

Template Parameters
proptoTrue if calculation is up to proportion (double-only terms dropped).
jacobian_adjust_transformTrue if the log absolute Jacobian determinant of inverse parameter transforms is added to the log probability.
MClass of model.
Parameters
[in]modelModel.
[in]params_rReal-valued parameters.
[in]params_iInteger-valued parameters.
[in,out]msgs

Definition at line 41 of file util.hpp.

template<bool jacobian_adjust_transform, class M >
double stan::model::log_prob_propto ( const M &  model,
Eigen::VectorXd &  params_r,
std::ostream *  msgs = 0 
)

Helper function to calculate log probability for double scalars up to a proportion.

This implementation wraps the double values in stan::agrad::var and calls the model's log_prob() function with propto=true and the specified parameter for applying the Jacobian adjustment for transformed parameters.

Template Parameters
proptoTrue if calculation is up to proportion (double-only terms dropped).
jacobian_adjust_transformTrue if the log absolute Jacobian determinant of inverse parameter transforms is added to the log probability.
MClass of model.
Parameters
[in]modelModel.
[in]params_rReal-valued parameters.
[in,out]msgs

Definition at line 133 of file util.hpp.

template<bool propto, bool jacobian_adjust_transform, class M >
int stan::model::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.

This shouldn't be necessary when using autodiff, but is useful for finding bugs in hand-written code (or agrad).

Template Parameters
proptoTrue if calculation is up to proportion (double-only terms dropped).
jacobian_adjust_transformTrue if the log absolute Jacobian determinant of inverse parameter transforms is added to the log probability.
MClass of model.
Parameters
modelModel.
params_rReal-valued parameter vector.
params_iInteger-valued parameter vector.
epsilonReal-valued scalar saying how much to perturb. Defaults to 1e-6.
errorReal-valued scalar saying how much error to allow. Defaults to 1e-6.
oOutput stream for messages. Defaults to std::cout.
msgsStream to which Stan programs write. Defaults to 0.
Returns
number of failed gradient comparisons versus allowed error, so 0 if all gradients pass

Definition at line 268 of file util.hpp.


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