Stan Math Library  2.12.0
reverse mode automatic differentiation
ibeta.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_REV_SCAL_FUN_IBETA_HPP
2 #define STAN_MATH_REV_SCAL_FUN_IBETA_HPP
3 
4 #include <stan/math/rev/core.hpp>
6 #include <boost/math/special_functions/digamma.hpp>
7 #include <boost/math/special_functions/gamma.hpp>
8 
9 namespace stan {
10  namespace math {
11 
12  namespace {
18  double ibeta_hypergeometric_helper(double a, double b, double z,
19  double precision = 1e-8,
20  double max_steps = 1000) {
21  double val = 0;
22  double diff = 1;
23  double k = 0;
24  double a_2 = a * a;
25  double bprod = 1;
26  while (std::abs(diff) > precision && ++k < max_steps) {
27  val += diff;
28  bprod *= b + k - 1.0;
29  diff = a_2 * std::pow(a + k, -2) * bprod * std::pow(z, k)
30  / boost::math::tgamma(k + 1);
31  }
32  return val;
33  }
34 
35  class ibeta_vvv_vari : public op_vvv_vari {
36  public:
37  ibeta_vvv_vari(vari* avi, vari* bvi, vari* xvi) :
38  op_vvv_vari(ibeta(avi->val_, bvi->val_, xvi->val_),
39  avi, bvi, xvi) {
40  }
41  void chain() {
42  double a = avi_->val_;
43  double b = bvi_->val_;
44  double c = cvi_->val_;
45 
46  using std::sin;
47  using std::pow;
48  using std::log;
50  using boost::math::tgamma;
52  using boost::math::ibeta;
53  avi_->adj_ += adj_ *
54  (log(c) - digamma(a) + digamma(a+b)) * val_
55  - tgamma(a) * tgamma(a+b) / tgamma(b) * pow(c, a)
56  / tgamma(1+a) / tgamma(1+a)
57  * ibeta_hypergeometric_helper(a, 1-b, c);
58  bvi_->adj_ += adj_ *
59  (tgamma(b) * tgamma(a+b) / tgamma(a) * pow(1-c, b)
60  * ibeta_hypergeometric_helper(b, 1-a, 1-c)
61  / tgamma(b+1) / tgamma(b+1)
62  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
63  cvi_->adj_ += adj_ *
64  boost::math::ibeta_derivative(a, b, c);
65  }
66  };
67  class ibeta_vvd_vari : public op_vvd_vari {
68  public:
69  ibeta_vvd_vari(vari* avi, vari* bvi, double x) :
70  op_vvd_vari(ibeta(avi->val_, bvi->val_, x), avi, bvi, x) {
71  }
72  void chain() {
73  double a = avi_->val_;
74  double b = bvi_->val_;
75  double c = cd_;
76 
77  using std::sin;
78  using std::pow;
79  using std::log;
81  using boost::math::tgamma;
83  using boost::math::ibeta;
84  avi_->adj_ += adj_ *
85  (log(c) - digamma(a) + digamma(a+b)) * val_ -
86  tgamma(a) * tgamma(a+b) / tgamma(b) * pow(c, a)
87  / tgamma(1+a) / tgamma(1+a)
88  * ibeta_hypergeometric_helper(a, 1-b, c);
89  bvi_->adj_ += adj_ *
90  (tgamma(b) * tgamma(a+b) / tgamma(a) * pow(1-c, b)
91  * ibeta_hypergeometric_helper(b, 1-a, 1-c)
92  / tgamma(b+1) / tgamma(b+1)
93  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
94  }
95  };
96  class ibeta_vdv_vari : public op_vdv_vari {
97  public:
98  ibeta_vdv_vari(vari* avi, double b, vari* xvi) :
99  op_vdv_vari(ibeta(avi->val_, b, xvi->val_), avi, b, xvi) {
100  }
101  void chain() {
102  double a = avi_->val_;
103  double b = bd_;
104  double c = cvi_->val_;
105 
106  using std::sin;
107  using std::pow;
108  using std::log;
110  using boost::math::tgamma;
111  using boost::math::digamma;
112  using boost::math::ibeta;
113  avi_->adj_ += adj_ *
114  (log(c) - digamma(a) + digamma(a+b)) * val_
115  - tgamma(a) * tgamma(a+b) / tgamma(b) * pow(c, a)
116  / tgamma(1+a) / tgamma(1+a)
117  * ibeta_hypergeometric_helper(a, 1-b, c);
118  cvi_->adj_ += adj_ *
119  boost::math::ibeta_derivative(a, b, c);
120  }
121  };
122  class ibeta_vdd_vari : public op_vdd_vari {
123  public:
124  ibeta_vdd_vari(vari* avi, double b, double x) :
125  op_vdd_vari(ibeta(avi->val_, b, x), avi, b, x) {
126  }
127  void chain() {
128  double a = avi_->val_;
129  double b = bd_;
130  double c = cd_;
131 
132  using std::sin;
133  using std::pow;
134  using std::log;
136  using boost::math::tgamma;
137  using boost::math::digamma;
138  using boost::math::ibeta;
139  avi_->adj_ += adj_ *
140  (log(c) - digamma(a) + digamma(a+b)) * val_
141  - tgamma(a) * tgamma(a+b) / tgamma(b) * pow(c, a)
142  / tgamma(1+a) / tgamma(1+a)
143  * ibeta_hypergeometric_helper(a, 1-b, c);
144  }
145  };
146  class ibeta_dvv_vari : public op_dvv_vari {
147  public:
148  ibeta_dvv_vari(double a, vari* bvi, vari* xvi) :
149  op_dvv_vari(ibeta(a, bvi->val_, xvi->val_), a, bvi, xvi) {
150  }
151  void chain() {
152  double a = ad_;
153  double b = bvi_->val_;
154  double c = cvi_->val_;
155 
156  using std::sin;
157  using std::pow;
158  using std::log;
160  using boost::math::tgamma;
161  using boost::math::digamma;
162  using boost::math::ibeta;
163  bvi_->adj_ += adj_ *
164  (tgamma(b) * tgamma(a+b) / tgamma(a) * pow(1-c, b)
165  * ibeta_hypergeometric_helper(b, 1-a, 1-c)
166  / tgamma(b+1) / tgamma(b+1)
167  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
168  cvi_->adj_ += adj_ *
169  boost::math::ibeta_derivative(a, b, c);
170  }
171  };
172  class ibeta_dvd_vari : public op_dvd_vari {
173  public:
174  ibeta_dvd_vari(double a, vari* bvi, double x) :
175  op_dvd_vari(ibeta(a, bvi->val_, x), a, bvi, x) {
176  }
177  void chain() {
178  double a = ad_;
179  double b = bvi_->val_;
180  double c = cd_;
181 
182  using std::sin;
183  using std::pow;
184  using std::log;
186  using boost::math::tgamma;
187  using boost::math::digamma;
188  using boost::math::ibeta;
189  bvi_->adj_ += adj_ *
190  (tgamma(b) * tgamma(a+b) / tgamma(a) * pow(1-c, b)
191  * ibeta_hypergeometric_helper(b, 1-a, 1-c)
192  / tgamma(b+1) / tgamma(b+1)
193  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
194  }
195  };
196  class ibeta_ddv_vari : public op_ddv_vari {
197  public:
198  ibeta_ddv_vari(double a, double b, vari* xvi) :
199  op_ddv_vari(ibeta(a, b, xvi->val_), a, b, xvi) {
200  }
201  void chain() {
202  double a = ad_;
203  double b = bd_;
204  double c = cvi_->val_;
205 
206  cvi_->adj_ += adj_ *
207  boost::math::ibeta_derivative(a, b, c);
208  }
209  };
210  }
211 
230  inline var ibeta(const var& a,
231  const var& b,
232  const var& x) {
233  return var(new ibeta_vvv_vari(a.vi_, b.vi_, x.vi_));
234  }
235 
236  }
237 }
238 #endif
fvar< T > abs(const fvar< T > &x)
Definition: abs.hpp:15
double ibeta(const double a, const double b, const double x)
The normalized incomplete beta function of a, b, and x.
Definition: ibeta.hpp:23
fvar< T > log(const fvar< T > &x)
Definition: log.hpp:14
Independent (input) and dependent (output) variables for gradients.
Definition: var.hpp:30
fvar< T > sin(const fvar< T > &x)
Definition: sin.hpp:12
vari * vi_
Pointer to the implementation of this variable.
Definition: var.hpp:42
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:94
var ibeta(const var &a, const var &b, const var &x)
The normalized incomplete beta function of a, b, and x.
Definition: ibeta.hpp:230
double pi()
Return the value of pi.
Definition: constants.hpp:85
double adj_
The adjoint of this variable, which is the partial derivative of this variable with respect to the ro...
Definition: vari.hpp:44
fvar< T > pow(const fvar< T > &x1, const fvar< T > &x2)
Definition: pow.hpp:17
fvar< T > tgamma(const fvar< T > &x)
Definition: tgamma.hpp:14
fvar< T > digamma(const fvar< T > &x)
Definition: digamma.hpp:15

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