Stan  1.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerator Friends Macros Pages
prob_grad_ad.hpp
Go to the documentation of this file.
1 #ifndef __STAN__MODEL__PROB_GRAD_AD_HPP__
2 #define __STAN__MODEL__PROB_GRAD_AD_HPP__
3 
4 #include <cstddef>
5 #include <utility>
6 #include <vector>
7 
8 #include <stan/agrad/agrad.hpp>
10 
11 namespace stan {
12 
13  namespace model {
14 
15  class prob_grad_ad : public prob_grad {
16  public:
17 
19  : prob_grad::prob_grad(num_params_r) {
20  }
21 
23  std::vector<std::pair<int,int> >& param_ranges_i)
24  : prob_grad::prob_grad(num_params_r,
25  param_ranges_i) {
26  }
27 
28  virtual ~prob_grad_ad() {
29  }
30 
31  virtual agrad::var log_prob(std::vector<agrad::var>& params_r,
32  std::vector<int>& params_i,
33  std::ostream* output_stream = 0) = 0;
34 
35  virtual double grad_log_prob(std::vector<double>& params_r,
36  std::vector<int>& params_i,
37  std::vector<double>& gradient,
38  std::ostream* output_stream = 0) {
39  std::vector<agrad::var> ad_params_r(num_params_r());
40  for (size_t i = 0; i < num_params_r(); ++i) {
41  agrad::var var_i(params_r[i]);
42  ad_params_r[i] = var_i;
43  }
44  agrad::var adLogProb = log_prob(ad_params_r,params_i,output_stream);
45  double val = adLogProb.val();
46  adLogProb.grad(ad_params_r,gradient);
47  return val;
48  }
49 
50  virtual double log_prob(std::vector<double>& params_r,
51  std::vector<int>& params_i,
52  std::ostream* output_stream = 0) {
53  std::vector<agrad::var> ad_params_r;
54  for (size_t i = 0; i < num_params_r(); ++i) {
55  agrad::var var_i(params_r[i]);
56  ad_params_r.push_back(var_i);
57  }
58  agrad::var adLogProb = log_prob(ad_params_r,params_i,output_stream);
59  double val = adLogProb.val();
61  return val;
62  }
63 
64  };
65 
66  }
67 }
68 
69 #endif
70 

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