1 #ifndef __STAN__MODEL__PROB_GRAD_HPP__
2 #define __STAN__MODEL__PROB_GRAD_HPP__
37 std::vector<std::pair<int,int> >& param_ranges_i)
80 virtual void init(std::vector<double>& params_r,
81 std::vector<int>& params_i) {
89 std::vector<int>& params_i,
90 std::vector<double>& gradient,
91 std::ostream* output_stream = 0) = 0;
93 virtual double log_prob(std::vector<double>& params_r,
94 std::vector<int>& params_i,
95 std::ostream* output_stream = 0) = 0;
113 std::vector<int>& params_i,
114 std::vector<double>& gradient,
115 std::vector<double>& hessian,
116 std::ostream* output_stream = 0) {
120 const double coefficients[order] = {1.0/12.0,-2.0/3.0,2.0/3.0,-1.0/12.0};
125 hessian.assign(params_r.size() * params_r.size(), 0);
126 std::vector<double> temp_grad(params_r.size());
127 std::vector<double> perturbed_params(params_r.begin(), params_r.end());
128 for (
size_t d = 0; d < params_r.size(); d++) {
129 double*
row = &hessian[d*params_r.size()];
130 for (
int i = 0; i < order; i++) {
131 perturbed_params[d] = params_r[d] + perturbations[i];
133 for (
size_t dd = 0; dd < params_r.size(); dd++) {
134 row[dd] += 0.5 * coefficients[i] * temp_grad[dd] /
epsilon;
135 hessian[d + dd*params_r.size()] += 0.5 * coefficients[i] * temp_grad[dd] /
epsilon;
138 perturbed_params[d] = params_r[d];
146 std::vector<double>& params_r,
147 std::vector<int>& params_i,
148 std::ostream* output_stream = 0) {
150 throw std::runtime_error (
"idx >= num_params_i()");
152 throw std::runtime_error (
"val <= param_range_i(idx) lower");
154 throw std::runtime_error (
"val >= param_range_i(idx) upper");
156 int original_val = params_i[idx];
158 double result =
log_prob(params_r,params_i,output_stream);
159 params_i[idx] = original_val;
178 std::vector<int>& params_i,
180 std::ostream* output_stream = 0) {
182 for (
size_t i = 0; i < params_i.size(); ++i)
183 writer.
write(params_i[i]);
184 for (
size_t i = 0; i < params_r.size(); ++i)
185 writer.
write(params_r[i]);
202 std::vector<int>& params_i,
203 std::vector<double>&
grad,
205 std::ostream* output_stream = 0) {
206 std::vector<double> perturbed(params_r);
208 grad.resize(params_r.size());
209 for (
size_t k = 0; k < params_r.size(); k++) {
211 double logp_plus =
log_prob(perturbed,params_i,output_stream);
212 perturbed[k] = params_r[k] -
epsilon;
213 double logp_minus =
log_prob(perturbed,params_i,output_stream);
214 double gradest = (logp_plus - logp_minus) / (2*
epsilon);
216 perturbed[k] = params_r[k];
237 std::vector<int>& params_i,
240 std::ostream& o = std::cout,
241 std::ostream* output_stream = 0) {
242 std::vector<double>
grad;
245 std::vector<double> grad_fd;
251 <<
" Log probability=" << lp
255 << std::setw(10) <<
"param idx"
256 << std::setw(16) <<
"value"
257 << std::setw(16) <<
"model"
258 << std::setw(16) <<
"finite diff"
259 << std::setw(16) <<
"error"
261 for (
size_t k = 0; k < params_r.size(); k++) {
262 o << std::setw(10) << k
263 << std::setw(16) << params_r[k]
264 << std::setw(16) <<
grad[k]
265 << std::setw(16) << grad_fd[k]
266 << std::setw(16) << (
grad[k] - grad_fd[k])
Writes Stan variables in comma-separated-value format to an output stream.
void newline()
Write a newline character.
void write(int n)
Write a value.
The prob_grad class represents densities with fixed numbers of discrete and scalar parameters and the...
std::vector< std::pair< int, int > > param_ranges_i__
virtual void init(std::vector< double > ¶ms_r, std::vector< int > ¶ms_i)
virtual size_t num_params_r()
virtual void write_csv(std::vector< double > ¶ms_r, std::vector< int > ¶ms_i, std::ostream &o, std::ostream *output_stream=0)
Write the parameters on a single line in CSV format.
std::pair< int, int > param_range_i(size_t idx)
virtual size_t num_params_i()
virtual double log_prob_star(size_t idx, int val, std::vector< double > ¶ms_r, std::vector< int > ¶ms_i, std::ostream *output_stream=0)
void set_num_params_r(size_t num_params_r)
prob_grad(size_t num_params_r, std::vector< std::pair< int, int > > ¶m_ranges_i)
prob_grad(size_t num_params_r)
void set_param_range_i_lower(size_t idx, int low)
int test_gradients(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 *output_stream=0)
Test the grad_log_prob() function's ability to produce accurate gradients using finite differences.
virtual double grad_log_prob(std::vector< double > ¶ms_r, std::vector< int > ¶ms_i, std::vector< double > &gradient, std::ostream *output_stream=0)=0
void setparam_ranges_i__(std::vector< std::pair< int, int > > param_ranges_i)
void finite_diff_grad(std::vector< double > ¶ms_r, std::vector< int > ¶ms_i, std::vector< double > &grad, double epsilon=1e-6, std::ostream *output_stream=0)
Compute the gradient using finite differences for the specified parameters, writing the result into t...
virtual double log_prob(std::vector< double > ¶ms_r, std::vector< int > ¶ms_i, std::ostream *output_stream=0)=0
virtual double grad_hess_log_prob(std::vector< double > ¶ms_r, std::vector< int > ¶ms_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.
void set_param_range_i_upper(size_t idx, int up)
int param_range_i_lower(size_t idx)
int param_range_i_upper(size_t idx)
var fabs(const var &a)
Return the absolute value of the variable (cmath).
static void grad(chainable *vi)
Compute the gradient for all variables starting from the specified root variable implementation.
double epsilon()
Return minimum positive number representable.
double e()
Return the base of the natural logarithm.
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.
Probability, optimization and sampling library.
Template specification of functions in std for Stan.