Stan Math Library  2.8.0
reverse mode automatic differentiation
 All Classes Namespaces Files Functions Variables Typedefs Enumerator Friends Macros Groups
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 <boost/math/special_functions/digamma.hpp>
5 #include <boost/math/special_functions/gamma.hpp>
6 #include <stan/math/rev/core.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
27  && ++k < max_steps
28  && !std::isnan(diff)) {
29  val += diff;
30  bprod *= b+k-1.0;
31  diff = a_2 * std::pow(a+k, -2) * bprod * std::pow(z, k)
32  / boost::math::tgamma(k+1);
33  }
34  return val;
35  }
36 
37  class ibeta_vvv_vari : public op_vvv_vari {
38  public:
39  ibeta_vvv_vari(vari* avi, vari* bvi, vari* xvi) :
40  op_vvv_vari(stan::math::ibeta(avi->val_, bvi->val_, xvi->val_),
41  avi, bvi, xvi) {
42  }
43  void chain() {
44  double a = avi_->val_;
45  double b = bvi_->val_;
46  double c = cvi_->val_;
47 
48  using std::sin;
49  using std::pow;
50  using std::log;
52  using boost::math::tgamma;
54  using boost::math::ibeta;
55  using stan::math::ibeta_hypergeometric_helper;
56  avi_->adj_ += adj_ *
57  (log(c) - digamma(a) + digamma(a+b)) * val_
58  - tgamma(a) * tgamma(a+b) / tgamma(b) * pow(c, a)
59  / tgamma(1+a) / tgamma(1+a)
60  * ibeta_hypergeometric_helper(a, 1-b, c);
61  bvi_->adj_ += adj_ *
62  (tgamma(b) * tgamma(a+b) / tgamma(a) * pow(1-c, b)
63  * ibeta_hypergeometric_helper(b, 1-a, 1-c)
64  / tgamma(b+1) / tgamma(b+1)
65  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
66  cvi_->adj_ += adj_ *
67  boost::math::ibeta_derivative(a, b, c);
68  }
69  };
70  class ibeta_vvd_vari : public op_vvd_vari {
71  public:
72  ibeta_vvd_vari(vari* avi, vari* bvi, double x) :
73  op_vvd_vari(stan::math::ibeta(avi->val_, bvi->val_, x), avi, bvi, x) {
74  }
75  void chain() {
76  double a = avi_->val_;
77  double b = bvi_->val_;
78  double c = cd_;
79 
80  using std::sin;
81  using std::pow;
82  using std::log;
84  using boost::math::tgamma;
86  using boost::math::ibeta;
87  using stan::math::ibeta_hypergeometric_helper;
88  avi_->adj_ += adj_ *
89  (log(c) - digamma(a) + digamma(a+b)) * val_ -
90  tgamma(a) * tgamma(a+b) / tgamma(b) * pow(c, a)
91  / tgamma(1+a) / tgamma(1+a)
92  * ibeta_hypergeometric_helper(a, 1-b, c);
93  bvi_->adj_ += adj_ *
94  (tgamma(b) * tgamma(a+b) / tgamma(a) * pow(1-c, b)
95  * ibeta_hypergeometric_helper(b, 1-a, 1-c)
96  / tgamma(b+1) / tgamma(b+1)
97  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
98  }
99  };
100  class ibeta_vdv_vari : public op_vdv_vari {
101  public:
102  ibeta_vdv_vari(vari* avi, double b, vari* xvi) :
103  op_vdv_vari(stan::math::ibeta(avi->val_, b, xvi->val_), avi, b, xvi) {
104  }
105  void chain() {
106  double a = avi_->val_;
107  double b = bd_;
108  double c = cvi_->val_;
109 
110  using std::sin;
111  using std::pow;
112  using std::log;
114  using boost::math::tgamma;
115  using boost::math::digamma;
116  using boost::math::ibeta;
117  using stan::math::ibeta_hypergeometric_helper;
118  avi_->adj_ += adj_ *
119  (log(c) - digamma(a) + digamma(a+b)) * val_
120  - tgamma(a) * tgamma(a+b) / tgamma(b) * pow(c, a)
121  / tgamma(1+a) / tgamma(1+a)
122  * ibeta_hypergeometric_helper(a, 1-b, c);
123  cvi_->adj_ += adj_ *
124  boost::math::ibeta_derivative(a, b, c);
125  }
126  };
127  class ibeta_vdd_vari : public op_vdd_vari {
128  public:
129  ibeta_vdd_vari(vari* avi, double b, double x) :
130  op_vdd_vari(stan::math::ibeta(avi->val_, b, x), avi, b, x) {
131  }
132  void chain() {
133  double a = avi_->val_;
134  double b = bd_;
135  double c = cd_;
136 
137  using std::sin;
138  using std::pow;
139  using std::log;
141  using boost::math::tgamma;
142  using boost::math::digamma;
143  using boost::math::ibeta;
144  using stan::math::ibeta_hypergeometric_helper;
145  avi_->adj_ += adj_ *
146  (log(c) - digamma(a) + digamma(a+b)) * val_
147  - tgamma(a) * tgamma(a+b) / tgamma(b) * pow(c, a)
148  / tgamma(1+a) / tgamma(1+a)
149  * ibeta_hypergeometric_helper(a, 1-b, c);
150  }
151  };
152  class ibeta_dvv_vari : public op_dvv_vari {
153  public:
154  ibeta_dvv_vari(double a, vari* bvi, vari* xvi) :
155  op_dvv_vari(stan::math::ibeta(a, bvi->val_, xvi->val_), a, bvi, xvi) {
156  }
157  void chain() {
158  double a = ad_;
159  double b = bvi_->val_;
160  double c = cvi_->val_;
161 
162  using std::sin;
163  using std::pow;
164  using std::log;
166  using boost::math::tgamma;
167  using boost::math::digamma;
168  using boost::math::ibeta;
169  using stan::math::ibeta_hypergeometric_helper;
170  bvi_->adj_ += adj_ *
171  (tgamma(b) * tgamma(a+b) / tgamma(a) * pow(1-c, b)
172  * ibeta_hypergeometric_helper(b, 1-a, 1-c)
173  / tgamma(b+1) / tgamma(b+1)
174  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
175  cvi_->adj_ += adj_ *
176  boost::math::ibeta_derivative(a, b, c);
177  }
178  };
179  class ibeta_dvd_vari : public op_dvd_vari {
180  public:
181  ibeta_dvd_vari(double a, vari* bvi, double x) :
182  op_dvd_vari(stan::math::ibeta(a, bvi->val_, x), a, bvi, x) {
183  }
184  void chain() {
185  double a = ad_;
186  double b = bvi_->val_;
187  double c = cd_;
188 
189  using std::sin;
190  using std::pow;
191  using std::log;
193  using boost::math::tgamma;
194  using boost::math::digamma;
195  using boost::math::ibeta;
196  using stan::math::ibeta_hypergeometric_helper;
197  bvi_->adj_ += adj_ *
198  (tgamma(b) * tgamma(a+b) / tgamma(a) * pow(1-c, b)
199  * ibeta_hypergeometric_helper(b, 1-a, 1-c)
200  / tgamma(b+1) / tgamma(b+1)
201  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
202  }
203  };
204  class ibeta_ddv_vari : public op_ddv_vari {
205  public:
206  ibeta_ddv_vari(double a, double b, vari* xvi) :
207  op_ddv_vari(stan::math::ibeta(a, b, xvi->val_), a, b, xvi) {
208  }
209  void chain() {
210  double a = ad_;
211  double b = bd_;
212  double c = cvi_->val_;
213 
214  cvi_->adj_ += adj_ *
215  boost::math::ibeta_derivative(a, b, c);
216  }
217  };
218  }
219 
238  inline var ibeta(const var& a,
239  const var& b,
240  const var& x) {
241  return var(new ibeta_vvv_vari(a.vi_, b.vi_, x.vi_));
242  }
243 
244  }
245 }
246 #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:15
Independent (input) and dependent (output) variables for gradients.
Definition: var.hpp:32
int isnan(const stan::math::var &a)
Checks if the given number is NaN.
Definition: std_isnan.hpp:18
fvar< T > sin(const fvar< T > &x)
Definition: sin.hpp:14
vari * vi_
Pointer to the implementation of this variable.
Definition: var.hpp:44
double e()
Return the base of the natural logarithm.
Definition: constants.hpp:95
var ibeta(const var &a, const var &b, const var &x)
The normalized incomplete beta function of a, b, and x.
Definition: ibeta.hpp:238
double pi()
Return the value of pi.
Definition: constants.hpp:86
double adj_
The adjoint of this variable, which is the partial derivative of this variable with respect to the ro...
Definition: vari.hpp:42
fvar< T > pow(const fvar< T > &x1, const fvar< T > &x2)
Definition: pow.hpp:18
fvar< T > tgamma(const fvar< T > &x)
Definition: tgamma.hpp:15
fvar< T > digamma(const fvar< T > &x)
Definition: digamma.hpp:16

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