Stan  1.0
probability, sampling & optimization
special_functions.hpp
Go to the documentation of this file.
1 #ifndef __STAN__AGRAD__AGRAD_SPECIAL_FUNCTIONS_HPP__
2 #define __STAN__AGRAD__AGRAD_SPECIAL_FUNCTIONS_HPP__
3 
6 
7 #include <boost/math/special_functions/acosh.hpp>
8 #include <boost/math/special_functions/asinh.hpp>
9 #include <boost/math/special_functions/atanh.hpp>
10 #include <boost/math/special_functions/digamma.hpp>
11 #include <boost/math/special_functions/hypot.hpp>
12 
13 namespace stan {
14 
15  namespace agrad {
16 
17  namespace {
18 
19  class lgamma_vari : public op_v_vari {
20  public:
21  lgamma_vari(double value, vari* avi) :
22  op_v_vari(value, avi) {
23  }
24  void chain() {
25  avi_->adj_ += adj_ * boost::math::digamma(avi_->val_);
26  }
27  };
28 
29  class tgamma_vari : public op_v_vari {
30  public:
31  tgamma_vari(vari* avi) :
32  op_v_vari(boost::math::tgamma(avi->val_), avi) {
33  }
34  void chain() {
35  avi_->adj_ += adj_ * val_ * boost::math::digamma(avi_->val_);
36  }
37  };
38 
39  class log1p_vari : public op_v_vari {
40  public:
41  log1p_vari(vari* avi) :
42  op_v_vari(stan::math::log1p(avi->val_),avi) {
43  }
44  void chain() {
45  avi_->adj_ += adj_ / (1 + avi_->val_);
46  }
47  };
48 
49  class log1m_vari : public op_v_vari {
50  public:
51  log1m_vari(vari* avi) :
52  op_v_vari(stan::math::log1p(-avi->val_),avi) {
53  }
54  void chain() {
55  avi_->adj_ += adj_ / (avi_->val_ - 1);
56  }
57  };
58 
59  class binary_log_loss_1_vari : public op_v_vari {
60  public:
61  binary_log_loss_1_vari(vari* avi) :
62  op_v_vari(-std::log(avi->val_),avi) {
63  }
64  void chain() {
65  avi_->adj_ -= adj_ / avi_->val_;
66  }
67  };
68 
69  class binary_log_loss_0_vari : public op_v_vari {
70  public:
71  binary_log_loss_0_vari(vari* avi) :
72  op_v_vari(-stan::math::log1p(-avi->val_),avi) {
73  }
74  void chain() {
75  avi_->adj_ += adj_ / (1.0 - avi_->val_);
76  }
77  };
78 
79  class fdim_vv_vari : public op_vv_vari {
80  public:
81  fdim_vv_vari(vari* avi, vari* bvi) :
82  op_vv_vari(avi->val_ - bvi->val_, avi, bvi) {
83  }
84  void chain() {
85  avi_->adj_ += adj_;
86  bvi_->adj_ -= adj_;
87  }
88  };
89 
90  class fdim_vd_vari : public op_v_vari {
91  public:
92  fdim_vd_vari(vari* avi, double b) :
93  op_v_vari(avi->val_ - b, avi) {
94  }
95  void chain() {
96  avi_->adj_ += adj_;
97  }
98  };
99 
100  class fdim_dv_vari : public op_v_vari {
101  public:
102  fdim_dv_vari(double a, vari* bvi) :
103  op_v_vari(a - bvi->val_, bvi) {
104  }
105  void chain() {
106  // avi_ is bvi argument to constructor
107  avi_->adj_ -= adj_;
108  }
109  };
110 
111  class fma_vvv_vari : public op_vvv_vari {
112  public:
113  fma_vvv_vari(vari* avi, vari* bvi, vari* cvi) :
114  op_vvv_vari(avi->val_ * bvi->val_ + cvi->val_,
115  avi,bvi,cvi) {
116  }
117  void chain() {
118  avi_->adj_ += adj_ * bvi_->val_;
119  bvi_->adj_ += adj_ * avi_->val_;
120  cvi_->adj_ += adj_;
121  }
122  };
123 
124  class fma_vvd_vari : public op_vv_vari {
125  public:
126  fma_vvd_vari(vari* avi, vari* bvi, double c) :
127  op_vv_vari(avi->val_ * bvi->val_ + c,
128  avi,bvi) {
129  }
130  void chain() {
131  avi_->adj_ += adj_ * bvi_->val_;
132  bvi_->adj_ += adj_ * avi_->val_;
133  }
134  };
135 
136  class fma_vdv_vari : public op_vdv_vari {
137  public:
138  fma_vdv_vari(vari* avi, double b, vari* cvi) :
139  op_vdv_vari(avi->val_ * b + cvi->val_,
140  avi,b,cvi) {
141  }
142  void chain() {
143  avi_->adj_ += adj_ * bd_;
144  cvi_->adj_ += adj_;
145  }
146  };
147 
148  class fma_vdd_vari : public op_vd_vari {
149  public:
150  fma_vdd_vari(vari* avi, double b, double c) :
151  op_vd_vari(avi->val_ * b + c,
152  avi,b) {
153  }
154  void chain() {
155  avi_->adj_ += adj_ * bd_;
156  }
157  };
158 
159  class fma_ddv_vari : public op_v_vari {
160  public:
161  fma_ddv_vari(double a, double b, vari* cvi) :
162  op_v_vari(a * b + cvi->val_,
163  cvi) {
164  }
165  void chain() {
166  // avi_ is cvi from constructor
167  avi_->adj_ += adj_;
168  }
169  };
170 
171  class inv_logit_vari : public op_v_vari {
172  public:
173  inv_logit_vari(vari* avi) :
174  op_v_vari(math::inv_logit(avi->val_),avi) {
175  }
176  void chain() {
177  avi_->adj_ += adj_ * val_ * (1.0 - val_);
178  }
179  };
180 
181  class acosh_vari : public op_v_vari {
182  public:
183  acosh_vari(vari* avi) :
184  op_v_vari(boost::math::acosh(avi->val_),avi) {
185  }
186  void chain() {
187  avi_->adj_ += adj_ / std::sqrt(avi_->val_ * avi_->val_ - 1.0);
188  }
189  };
190 
191  class asinh_vari : public op_v_vari {
192  public:
193  asinh_vari(vari* avi) :
194  op_v_vari(boost::math::asinh(avi->val_),avi) {
195  }
196  void chain() {
197  avi_->adj_ += adj_ / std::sqrt(avi_->val_ * avi_->val_ + 1.0);
198  }
199  };
200 
201  class atanh_vari : public op_v_vari {
202  public:
203  atanh_vari(vari* avi) :
204  op_v_vari(boost::math::atanh(avi->val_),avi) {
205  }
206  void chain() {
207  avi_->adj_ += adj_ / (1.0 - avi_->val_ * avi_->val_);
208  }
209  };
210 
211  const double TWO_OVER_SQRT_PI = 2.0 / std::sqrt(boost::math::constants::pi<double>());
212 
213  class erf_vari : public op_v_vari {
214  public:
215  erf_vari(vari* avi) :
216  op_v_vari(boost::math::erf(avi->val_),avi) {
217  }
218  void chain() {
219  avi_->adj_ += adj_ * TWO_OVER_SQRT_PI * std::exp(- avi_->val_ * avi_->val_);
220  }
221  };
222 
223  const double NEG_TWO_OVER_SQRT_PI = - TWO_OVER_SQRT_PI;
224 
225  class erfc_vari : public op_v_vari {
226  public:
227  erfc_vari(vari* avi) :
228  op_v_vari(boost::math::erfc(avi->val_),avi) {
229  }
230  void chain() {
231  avi_->adj_ += adj_ * NEG_TWO_OVER_SQRT_PI * std::exp(- avi_->val_ * avi_->val_);
232  }
233  };
234 
235  const double LOG_2 = std::log(2.0);
236 
237  class exp2_vari : public op_v_vari {
238  public:
239  exp2_vari(vari* avi) :
240  op_v_vari(std::pow(2.0,avi->val_),avi) {
241  }
242  void chain() {
243  avi_->adj_ += adj_ * val_ * LOG_2;
244  }
245  };
246 
247  class expm1_vari : public op_v_vari {
248  public:
249  expm1_vari(vari* avi) :
250  op_v_vari(std::exp(avi->val_) - 1.0,avi) {
251  }
252  void chain() {
253  avi_->adj_ += adj_ * val_;
254  }
255  };
256 
257  class hypot_vv_vari : public op_vv_vari {
258  public:
259  hypot_vv_vari(vari* avi, vari* bvi) :
260  op_vv_vari(boost::math::hypot(avi->val_,bvi->val_),
261  avi,bvi) {
262  }
263  void chain() {
264  avi_->adj_ += adj_ * avi_->val_ / val_;
265  bvi_->adj_ += adj_ * bvi_->val_ / val_;
266  }
267  };
268 
269  class hypot_vd_vari : public op_v_vari {
270  public:
271  hypot_vd_vari(vari* avi, double b) :
272  op_v_vari(boost::math::hypot(avi->val_,b),
273  avi) {
274  }
275  void chain() {
276  avi_->adj_ += adj_ * avi_->val_ / val_;
277  }
278  };
279 
280  const double LOG2 = std::log(2.0);
281 
282  class log2_vari : public op_v_vari {
283  public:
284  log2_vari(vari* avi) :
285  op_v_vari(stan::math::log2(avi->val_),avi) {
286  }
287  void chain() {
288  avi_->adj_ += adj_ / (LOG2 * avi_->val_);
289  }
290  };
291 
292  class cbrt_vari : public op_v_vari {
293  public:
294  cbrt_vari(vari* avi) :
295  op_v_vari(boost::math::cbrt(avi->val_),avi) {
296  }
297  void chain() {
298  avi_->adj_ += adj_ / (3.0 * val_ * val_);
299  }
300  };
301 
302  // derivative 0 almost everywhere
303  class round_vari : public vari {
304  public:
305  round_vari(vari* avi) :
306  vari(boost::math::round(avi->val_)) {
307  }
308  };
309 
310  // derivative 0 almost everywhere
311  class trunc_vari : public vari {
312  public:
313  trunc_vari(vari* avi) :
314  vari(boost::math::trunc(avi->val_)) {
315  }
316  };
317 
318  class inv_cloglog_vari : public op_v_vari {
319  public:
320  inv_cloglog_vari(vari* avi) :
321  op_v_vari(stan::math::inv_cloglog(avi->val_), avi) {
322  }
323  void chain() {
324  avi_->adj_ -= adj_ * std::exp(avi_->val_ - std::exp(avi_->val_));
325  }
326  };
327 
328  class Phi_vari : public op_v_vari {
329  public:
330  Phi_vari(vari* avi) :
331  op_v_vari(stan::math::Phi(avi->val_), avi) {
332  }
333  void chain() {
334  static const double NEG_HALF = -0.5;
335  static const double INV_SQRT_TWO_PI
336  = 1.0 / std::sqrt(2.0 * boost::math::constants::pi<double>());
337  avi_->adj_ += adj_ * INV_SQRT_TWO_PI * std::exp(NEG_HALF * avi_->val_ * avi_->val_);
338  }
339  };
340 
341  inline double calculate_chain(const double& x, const double& val) {
342  return std::exp(x - val); // works out to inv_logit(x)
343  }
344 
345  double log_sum_exp_as_double(const std::vector<var>& x) {
346  using std::numeric_limits;
347  using std::exp;
348  using std::log;
349  double max = -numeric_limits<double>::infinity();
350  for (size_t i = 0; i < x.size(); ++i)
351  if (x[i] > max)
352  max = x[i].val();
353  double sum = 0.0;
354  for (size_t i = 0; i < x.size(); ++i)
355  if (x[i] != -numeric_limits<double>::infinity())
356  sum += exp(x[i].val() - max);
357  return max + log(sum);
358  }
359 
360  class log1p_exp_v_vari : public op_v_vari {
361  public:
362  log1p_exp_v_vari(vari* avi) :
363  op_v_vari(stan::math::log1p_exp(avi->val_),
364  avi) {
365  }
366  void chain() {
367  avi_->adj_ += adj_ * calculate_chain(avi_->val_, val_);
368  }
369  };
370  class log_sum_exp_vv_vari : public op_vv_vari {
371  public:
372  log_sum_exp_vv_vari(vari* avi, vari* bvi) :
373  op_vv_vari(stan::math::log_sum_exp(avi->val_, bvi->val_),
374  avi, bvi) {
375  }
376  void chain() {
377  avi_->adj_ += adj_ * calculate_chain(avi_->val_, val_);
378  bvi_->adj_ += adj_ * calculate_chain(bvi_->val_, val_);
379  }
380  };
381  class log_sum_exp_vd_vari : public op_vd_vari {
382  public:
383  log_sum_exp_vd_vari(vari* avi, double b) :
384  op_vd_vari(stan::math::log_sum_exp(avi->val_, b),
385  avi, b) {
386  }
387  void chain() {
388  avi_->adj_ += adj_ * calculate_chain(avi_->val_, val_);
389  }
390  };
391  class log_sum_exp_dv_vari : public op_dv_vari {
392  public:
393  log_sum_exp_dv_vari(double a, vari* bvi) :
394  op_dv_vari(stan::math::log_sum_exp(a, bvi->val_),
395  a, bvi) {
396  }
397  void chain() {
398  bvi_->adj_ += adj_ * calculate_chain(bvi_->val_, val_);
399  }
400  };
401  class log_sum_exp_vector_vari : public op_vector_vari {
402  public:
403  log_sum_exp_vector_vari(const std::vector<var>& x) :
404  op_vector_vari(log_sum_exp_as_double(x), x) {
405  }
406  void chain() {
407  for (size_t i = 0; i < size_; ++i) {
408  vis_[i]->adj_ += adj_ * calculate_chain(vis_[i]->val_, val_);
409  }
410  }
411  };
412 
413  class square_vari : public op_v_vari {
414  public:
415  square_vari(vari* avi) :
416  op_v_vari(avi->val_ * avi->val_,avi) {
417  }
418  void chain() {
419  avi_->adj_ += adj_ * 2.0 * avi_->val_;
420  }
421  };
422 
423  class multiply_log_vv_vari : public op_vv_vari {
424  public:
425  multiply_log_vv_vari(vari* avi, vari* bvi) :
426  op_vv_vari(stan::math::multiply_log(avi->val_,bvi->val_),avi,bvi) {
427  }
428  void chain() {
429  using std::log;
430  avi_->adj_ += adj_ * log(bvi_->val_);
431  if (bvi_->val_ == 0.0 && avi_->val_ == 0)
432  bvi_->adj_ += adj_ * std::numeric_limits<double>::infinity();
433  else
434  bvi_->adj_ += adj_ * avi_->val_ / bvi_->val_;
435  }
436  };
437  class multiply_log_vd_vari : public op_vd_vari {
438  public:
439  multiply_log_vd_vari(vari* avi, double b) :
440  op_vd_vari(stan::math::multiply_log(avi->val_,b),avi,b) {
441  }
442  void chain() {
443  using std::log;
444  avi_->adj_ += adj_ * log(bd_);
445  }
446  };
447  class multiply_log_dv_vari : public op_dv_vari {
448  public:
449  multiply_log_dv_vari(double a, vari* bvi) :
450  op_dv_vari(stan::math::multiply_log(a,bvi->val_),a,bvi) {
451  }
452  void chain() {
453  if (bvi_->val_ == 0.0 && ad_ == 0.0)
454  bvi_->adj_ += adj_ * std::numeric_limits<double>::infinity();
455  else
456  bvi_->adj_ += adj_ * ad_ / bvi_->val_;
457  }
458  };
459  namespace {
465  double ibeta_hypergeometric_helper(double a, double b, double z, double precision=1e-8, double max_steps=1000) {
466  double val=0;
467  double diff=1;
468  double k=0;
469  double a_2 = a*a;
470  double bprod = 1;
471  while (std::abs(diff) > precision && (++k < max_steps) && !std::isnan(diff)) {
472  val += diff;
473  bprod *= b+k-1.0;
474  diff = a_2*std::pow(a+k,-2)*bprod*std::pow(z,k)/boost::math::tgamma(k+1);
475  }
476  return val;
477  }
478  }
479  class ibeta_vvv_vari : public op_vvv_vari {
480  public:
481  ibeta_vvv_vari(vari* avi, vari* bvi, vari* xvi) :
482  op_vvv_vari(stan::math::ibeta(avi->val_,bvi->val_,xvi->val_),avi,bvi,xvi) {
483  }
484  void chain() {
485  double a = avi_->val_;
486  double b = bvi_->val_;
487  double c = cvi_->val_;
488 
489  using std::sin;
490  using std::pow;
491  using std::log;
493  using boost::math::tgamma;
494  using boost::math::digamma;
495  using boost::math::ibeta;
496  using stan::agrad::ibeta_hypergeometric_helper;
497  avi_->adj_ += adj_ *
498  (log(c) - digamma(a) + digamma(a+b)) * val_ -
499  tgamma(a)*tgamma(a+b)/tgamma(b) * pow(c,a) / tgamma(1+a) / tgamma(1+a) * ibeta_hypergeometric_helper(a, 1-b, c);
500  bvi_->adj_ += adj_ *
501  (tgamma(b)*tgamma(a+b)/tgamma(a)*pow(1-c,b) * ibeta_hypergeometric_helper(b,1-a,1-c)/tgamma(b+1)/tgamma(b+1)
502  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
503  cvi_->adj_ += adj_ *
504  boost::math::ibeta_derivative(a,b,c);
505  }
506  };
507  class ibeta_vvd_vari : public op_vvd_vari {
508  public:
509  ibeta_vvd_vari(vari* avi, vari* bvi, double x) :
510  op_vvd_vari(stan::math::ibeta(avi->val_,bvi->val_,x),avi,bvi,x) {
511  }
512  void chain() {
513  double a = avi_->val_;
514  double b = bvi_->val_;
515  double c = cd_;
516 
517  using std::sin;
518  using std::pow;
519  using std::log;
521  using boost::math::tgamma;
522  using boost::math::digamma;
523  using boost::math::ibeta;
524  using stan::agrad::ibeta_hypergeometric_helper;
525  avi_->adj_ += adj_ *
526  (log(c) - digamma(a) + digamma(a+b)) * val_ -
527  tgamma(a)*tgamma(a+b)/tgamma(b) * pow(c,a) / tgamma(1+a) / tgamma(1+a) * ibeta_hypergeometric_helper(a, 1-b, c);
528  bvi_->adj_ += adj_ *
529  (tgamma(b)*tgamma(a+b)/tgamma(a)*pow(1-c,b) * ibeta_hypergeometric_helper(b,1-a,1-c)/tgamma(b+1)/tgamma(b+1)
530  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
531  }
532  };
533  class ibeta_vdv_vari : public op_vdv_vari {
534  public:
535  ibeta_vdv_vari(vari* avi, double b, vari* xvi) :
536  op_vdv_vari(stan::math::ibeta(avi->val_,b,xvi->val_),avi,b,xvi) {
537  }
538  void chain() {
539  double a = avi_->val_;
540  double b = bd_;
541  double c = cvi_->val_;
542 
543  using std::sin;
544  using std::pow;
545  using std::log;
547  using boost::math::tgamma;
548  using boost::math::digamma;
549  using boost::math::ibeta;
550  using stan::agrad::ibeta_hypergeometric_helper;
551  avi_->adj_ += adj_ *
552  (log(c) - digamma(a) + digamma(a+b)) * val_ -
553  tgamma(a)*tgamma(a+b)/tgamma(b) * pow(c,a) / tgamma(1+a) / tgamma(1+a) * ibeta_hypergeometric_helper(a, 1-b, c);
554  cvi_->adj_ += adj_ *
555  boost::math::ibeta_derivative(a,b,c);
556  }
557  };
558  class ibeta_vdd_vari : public op_vdd_vari {
559  public:
560  ibeta_vdd_vari(vari* avi, double b, double x) :
561  op_vdd_vari(stan::math::ibeta(avi->val_,b,x),avi,b,x) {
562  }
563  void chain() {
564  double a = avi_->val_;
565  double b = bd_;
566  double c = cd_;
567 
568  using std::sin;
569  using std::pow;
570  using std::log;
572  using boost::math::tgamma;
573  using boost::math::digamma;
574  using boost::math::ibeta;
575  using stan::agrad::ibeta_hypergeometric_helper;
576  avi_->adj_ += adj_ *
577  (log(c) - digamma(a) + digamma(a+b)) * val_ -
578  tgamma(a)*tgamma(a+b)/tgamma(b) * pow(c,a) / tgamma(1+a) / tgamma(1+a) * ibeta_hypergeometric_helper(a, 1-b, c);
579  }
580  };
581  class ibeta_dvv_vari : public op_dvv_vari {
582  public:
583  ibeta_dvv_vari(double a, vari* bvi, vari* xvi) :
584  op_dvv_vari(stan::math::ibeta(a,bvi->val_,xvi->val_),a,bvi,xvi) {
585  }
586  void chain() {
587  double a = ad_;
588  double b = bvi_->val_;
589  double c = cvi_->val_;
590 
591  using std::sin;
592  using std::pow;
593  using std::log;
595  using boost::math::tgamma;
596  using boost::math::digamma;
597  using boost::math::ibeta;
598  using stan::agrad::ibeta_hypergeometric_helper;
599  bvi_->adj_ += adj_ *
600  (tgamma(b)*tgamma(a+b)/tgamma(a)*pow(1-c,b) * ibeta_hypergeometric_helper(b,1-a,1-c)/tgamma(b+1)/tgamma(b+1)
601  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
602  cvi_->adj_ += adj_ *
603  boost::math::ibeta_derivative(a,b,c);
604  }
605  };
606  class ibeta_dvd_vari : public op_dvd_vari {
607  public:
608  ibeta_dvd_vari(double a, vari* bvi, double x) :
609  op_dvd_vari(stan::math::ibeta(a,bvi->val_,x),a,bvi,x) {
610  }
611  void chain() {
612  double a = ad_;
613  double b = bvi_->val_;
614  double c = cd_;
615 
616  using std::sin;
617  using std::pow;
618  using std::log;
620  using boost::math::tgamma;
621  using boost::math::digamma;
622  using boost::math::ibeta;
623  using stan::agrad::ibeta_hypergeometric_helper;
624  bvi_->adj_ += adj_ *
625  (tgamma(b)*tgamma(a+b)/tgamma(a)*pow(1-c,b) * ibeta_hypergeometric_helper(b,1-a,1-c)/tgamma(b+1)/tgamma(b+1)
626  + ibeta(b, a, 1-c) * (digamma(b) - digamma(a+b) - log(1-c)));
627  }
628  };
629  class ibeta_ddv_vari : public op_ddv_vari {
630  public:
631  ibeta_ddv_vari(double a, double b, vari* xvi) :
632  op_ddv_vari(stan::math::ibeta(a,b,xvi->val_),a,b,xvi) {
633  }
634  void chain() {
635  double a = ad_;
636  double b = bd_;
637  double c = cvi_->val_;
638 
639  cvi_->adj_ += adj_ *
640  boost::math::ibeta_derivative(a,b,c);
641  }
642  };
643  }
644 
645 
648 
661  inline var acosh(const stan::agrad::var& a) {
662  return var(new acosh_vari(a.vi_));
663  }
664 
677  inline var asinh(const stan::agrad::var& a) {
678  return var(new asinh_vari(a.vi_));
679  }
680 
693  inline var atanh(const stan::agrad::var& a) {
694  return var(new atanh_vari(a.vi_));
695  }
696 
709  inline var erf(const stan::agrad::var& a) {
710  return var(new erf_vari(a.vi_));
711  }
712 
725  inline var erfc(const stan::agrad::var& a) {
726  return var(new erfc_vari(a.vi_));
727  }
728 
741  inline var exp2(const stan::agrad::var& a) {
742  return var(new exp2_vari(a.vi_));
743  }
744 
757  inline var expm1(const stan::agrad::var& a) {
758  return var(new expm1_vari(a.vi_));
759  }
760 
771  inline var lgamma(const stan::agrad::var& a) {
772  double lgamma_a = boost::math::lgamma(a.val());
773  return var(new lgamma_vari(lgamma_a, a.vi_));
774  }
775 
786  inline var log1p(const stan::agrad::var& a) {
787  return var(new log1p_vari(a.vi_));
788  }
789 
800  inline var log1m(const stan::agrad::var& a) {
801  return var(new log1m_vari(a.vi_));
802  }
803 
824  inline var fma(const stan::agrad::var& a,
825  const stan::agrad::var& b,
826  const stan::agrad::var& c) {
827  return var(new fma_vvv_vari(a.vi_,b.vi_,c.vi_));
828  }
829 
848  inline var fma(const stan::agrad::var& a,
849  const stan::agrad::var& b,
850  const double& c) {
851  return var(new fma_vvd_vari(a.vi_,b.vi_,c));
852  }
853 
872  inline var fma(const stan::agrad::var& a,
873  const double& b,
874  const stan::agrad::var& c) {
875  return var(new fma_vdv_vari(a.vi_,b,c.vi_));
876  }
877 
894  inline var fma(const stan::agrad::var& a,
895  const double& b,
896  const double& c) {
897  return var(new fma_vdd_vari(a.vi_,b,c));
898  }
899 
916  inline var fma(const double& a,
917  const stan::agrad::var& b,
918  const double& c) {
919  return var(new fma_vdd_vari(b.vi_,a,c));
920  }
921 
938  inline var fma(const double& a,
939  const double& b,
940  const stan::agrad::var& c) {
941  return var(new fma_ddv_vari(a,b,c.vi_));
942  }
943 
962  inline var fma(const double& a,
963  const stan::agrad::var& b,
964  const stan::agrad::var& c) {
965  return var(new fma_vdv_vari(b.vi_,a,c.vi_)); // a-b symmetry
966  }
967 
968 
986  inline var fmax(const stan::agrad::var& a,
987  const stan::agrad::var& b) {
988  return a.vi_->val_ > b.vi_->val_ ? a : b;
989  }
990 
1007  inline var fmax(const stan::agrad::var& a,
1008  const double& b) {
1009  return a.vi_->val_ >= b ? a : var(b);
1010  }
1011 
1028  inline var fmax(const double& a,
1029  const stan::agrad::var& b) {
1030  return a > b.vi_->val_ ? var(a) : b;
1031  }
1032 
1046  inline var fmin(const stan::agrad::var& a,
1047  const stan::agrad::var& b) {
1048  return a.vi_->val_ < b.vi_->val_ ? a : b;
1049  }
1050 
1065  inline var fmin(const stan::agrad::var& a,
1066  const double& b) {
1067  return a.vi_->val_ <= b ? a : var(b);
1068  }
1069 
1086  inline var fmin(const double& a,
1087  const stan::agrad::var& b) {
1088  return a < b.vi_->val_ ? var(a) : b;
1089  }
1090 
1107  inline var hypot(const stan::agrad::var& a,
1108  const stan::agrad::var& b) {
1109  return var(new hypot_vv_vari(a.vi_,b.vi_));
1110  }
1111 
1126  inline var hypot(const stan::agrad::var& a,
1127  const double& b) {
1128  return var(new hypot_vd_vari(a.vi_,b));
1129  }
1130 
1145  inline var hypot(const double& a,
1146  const stan::agrad::var& b) {
1147  return var(new hypot_vd_vari(b.vi_,a));
1148  }
1149 
1162  inline var log2(const stan::agrad::var& a) {
1163  return var(new log2_vari(a.vi_));
1164  }
1165 
1178  inline var cbrt(const stan::agrad::var& a) {
1179  return var(new cbrt_vari(a.vi_));
1180  }
1181 
1196  inline var round(const stan::agrad::var& a) {
1197  return var(new round_vari(a.vi_));
1198  }
1199 
1213  inline var trunc(const stan::agrad::var& a) {
1214  return var(new trunc_vari(a.vi_));
1215  }
1216 
1217 
1241  inline var fdim(const stan::agrad::var& a,
1242  const stan::agrad::var& b) {
1243  if (a.vi_->val_ > b.vi_->val_)
1244  return var(new fdim_vv_vari(a.vi_,b.vi_));
1245  else
1246  return var(new vari(0.0));
1247  }
1248 
1266  inline var fdim(const double& a,
1267  const stan::agrad::var& b) {
1268  return a > b.vi_->val_
1269  ? var(new fdim_dv_vari(a,b.vi_))
1270  : var(new vari(0.0));
1271  }
1272 
1289  inline var fdim(const stan::agrad::var& a,
1290  const double& b) {
1291  return a.vi_->val_ > b
1292  ? var(new fdim_vd_vari(a.vi_,b))
1293  : var(new vari(0.0));
1294  }
1295 
1312  inline var tgamma(const stan::agrad::var& a) {
1313  return var(new tgamma_vari(a.vi_));
1314  }
1315 
1332  inline var step(const stan::agrad::var& a) {
1333  return var(new vari(a.vi_->val_ < 0.0 ? 0.0 : 1.0));
1334  }
1335 
1350  inline var inv_cloglog(const stan::agrad::var& a) {
1351  return var(new inv_cloglog_vari(a.vi_));
1352  }
1353 
1366  inline var Phi(const stan::agrad::var& a) {
1367  return var(new Phi_vari(a.vi_));
1368  }
1369 
1382  inline var inv_logit(const stan::agrad::var& a) {
1383  return var(new inv_logit_vari(a.vi_));
1384  }
1385 
1401  inline var log_loss(const int& y,
1402  const stan::agrad::var& y_hat) {
1403  return y == 0
1404  ? var(new binary_log_loss_0_vari(y_hat.vi_))
1405  : var(new binary_log_loss_1_vari(y_hat.vi_));
1406  }
1407 
1412  inline var log1p_exp(const stan::agrad::var& a) {
1413  return var(new log1p_exp_v_vari(a.vi_));
1414  }
1415 
1420  const stan::agrad::var& b) {
1421  return var(new log_sum_exp_vv_vari(a.vi_, b.vi_));
1422  }
1427  const double& b) {
1428  return var(new log_sum_exp_vd_vari(a.vi_, b));
1429  }
1433  inline var log_sum_exp(const double& a,
1434  const stan::agrad::var& b) {
1435  return var(new log_sum_exp_dv_vari(a, b.vi_));
1436  }
1440  inline var log_sum_exp(const std::vector<var>& x) {
1441  return var(new log_sum_exp_vector_vari(x));
1442  }
1443 
1453  inline var square(const var& x) {
1454  return var(new square_vari(x.vi_));
1455  }
1456 
1457  // OTHER FUNCTIONS: stan/math/special_functions.hpp implementations
1470  inline var multiply_log(const var& a, const var& b) {
1471  return var(new multiply_log_vv_vari(a.vi_,b.vi_));
1472  }
1483  inline var multiply_log(const var& a, const double b) {
1484  return var(new multiply_log_vd_vari(a.vi_,b));
1485  }
1497  inline var multiply_log(const double a, const var& b) {
1498  if (a == 1.0)
1499  return log(b);
1500  return var(new multiply_log_dv_vari(a,b.vi_));
1501  }
1502 
1503 
1512  inline var if_else(bool c, const var& y_true, const var&y_false) {
1513  return c ? y_true : y_false;
1514  }
1524  inline var if_else(bool c, double y_true, const var& y_false) {
1525  if (c)
1526  return var(y_true);
1527  else
1528  return y_false;
1529  }
1539  inline var if_else(bool c, const var& y_true, const double y_false) {
1540  if (c)
1541  return y_true;
1542  else
1543  return var(y_false);
1544  }
1545 
1563  inline var ibeta(const var& a,
1564  const var& b,
1565  const var& x) {
1566  return var(new ibeta_vvv_vari(a.vi_, b.vi_, x.vi_));
1567  }
1568 
1582  inline double value_of(const agrad::var& v) {
1583  return v.vi_->val_;
1584  }
1585 
1592  inline int as_bool(const agrad::var& v) {
1593  return 0.0 != v.vi_->val_;
1594  }
1595 
1596  } // namespace math
1597 
1598 } // namespace stan
1599 
1600 
1601 #endif
double cd_
Definition: agrad.hpp:668
vari ** vis_
Definition: agrad.hpp:752
double bd_
Definition: agrad.hpp:629
double ad_
Definition: agrad.hpp:640
vari * avi_
Definition: agrad.hpp:606
vari * cvi_
Definition: agrad.hpp:654
const size_t size_
Definition: agrad.hpp:751
vari * bvi_
Definition: agrad.hpp:617
Independent (input) and dependent (output) variables for gradients.
Definition: agrad.hpp:214
double val() const
Return the value of this variable.
Definition: agrad.hpp:396
vari * vi_
Pointer to the implementation of this variable.
Definition: agrad.hpp:227
The variable implementation base class.
Definition: agrad.hpp:104
const double val_
The value of this variable.
Definition: agrad.hpp:113
Reimplementing boost functionality.
var inv_logit(const stan::agrad::var &a)
The inverse logit function for variables (stan).
var expm1(const stan::agrad::var &a)
The exponentiation of the specified variable minus 1 (C99).
int as_bool(const agrad::var &v)
Return 1 if the argument is unequal to zero and 0 otherwise.
var multiply_log(const var &a, const var &b)
Return the value of a*log(b).
var asinh(const stan::agrad::var &a)
The inverse hyperbolic sine function for variables (C99).
var log1p_exp(const stan::agrad::var &a)
Return the log of 1 plus the exponential of the specified variable.
var sqrt(const var &a)
Return the square root of the specified variable (cmath).
Definition: agrad.hpp:1761
var square(const var &x)
Return the square of the input variable.
var fma(const stan::agrad::var &a, const stan::agrad::var &b, const stan::agrad::var &c)
The fused multiply-add function for three variables (C99).
var if_else(bool c, const var &y_true, const var &y_false)
If the specified condition is true, return the first variable, otherwise return the second variable.
var abs(const var &a)
Return the absolute value of the variable (std).
Definition: agrad.hpp:2149
var exp2(const stan::agrad::var &a)
Exponentiation base 2 function for variables (C99).
var step(const stan::agrad::var &a)
Return the step, or heaviside, function applied to the specified variable (stan).
var log1m(const stan::agrad::var &a)
The log (1 - x) function for variables.
var ibeta(const var &a, const var &b, const var &x)
The normalized incomplete beta function of a, b, and x.
var acosh(const stan::agrad::var &a)
The inverse hyperbolic cosine function for variables (C99).
var cbrt(const stan::agrad::var &a)
Returns the cube root of the specified variable (C99).
var erfc(const stan::agrad::var &a)
The complementary error function for variables (C99).
var log1p(const stan::agrad::var &a)
The log (1 + x) function for variables (C99).
var fmax(const stan::agrad::var &a, const stan::agrad::var &b)
Returns the maximum of the two variable arguments (C99).
var log_loss(const int &y, const stan::agrad::var &y_hat)
The log loss function for variables (stan).
var sin(const var &a)
Return the sine of a radian-scaled variable (cmath).
Definition: agrad.hpp:1847
var fdim(const stan::agrad::var &a, const stan::agrad::var &b)
Return the positive difference between the first variable's the value and the second's (C99).
double value_of(const agrad::var &v)
Return the value of the specified variable.
var lgamma(const stan::agrad::var &a)
The log gamma function for variables (C99).
var erf(const stan::agrad::var &a)
The error function for variables (C99).
var Phi(const stan::agrad::var &a)
The unit normal cumulative density function for variables (stan).
var pow(const var &base, const var &exponent)
Return the base raised to the power of the exponent (cmath).
Definition: agrad.hpp:1778
var log_sum_exp(const stan::agrad::var &a, const stan::agrad::var &b)
Returns the log sum of exponentials.
var log(const var &a)
Return the natural log of the specified variable (cmath).
Definition: agrad.hpp:1730
var log2(const stan::agrad::var &a)
Returns the base 2 logarithm of the specified variable (C99).
var hypot(const stan::agrad::var &a, const stan::agrad::var &b)
Returns the length of the hypoteneuse of a right triangle with sides of the specified lengths (C99).
var tgamma(const stan::agrad::var &a)
Return the Gamma function applied to the specified variable (C99).
var round(const stan::agrad::var &a)
Returns the rounded form of the specified variable (C99).
var sum(const Eigen::Matrix< var, R, C > &m)
Returns the sum of the coefficients of the specified matrix, column vector or row vector.
Definition: matrix.hpp:837
var trunc(const stan::agrad::var &a)
Returns the truncatation of the specified variable (C99).
var atanh(const stan::agrad::var &a)
The inverse hyperbolic tangent function for variables (C99).
var exp(const var &a)
Return the exponentiation of the specified variable (cmath).
Definition: agrad.hpp:1716
var fmin(const stan::agrad::var &a, const stan::agrad::var &b)
Returns the minimum of the two variable arguments (C99).
var inv_cloglog(const stan::agrad::var &a)
Return the inverse complementary log-log function applied specified variable (stan).
double e()
Return the base of the natural logarithm.
double ibeta(const double &a, const double &b, const double &x)
The normalized incomplete beta function of a, b, and x.
const double LOG_2
The natural logarithm of 2, .
Definition: constants.hpp:26
bool check_greater_or_equal(const char *function, const T_y &y, const T_low &low, const char *name, T_result *result, const Policy &)
int max(const std::vector< int > &x)
Returns the maximum coefficient in the specified column vector.
Definition: matrix.hpp:968
bool check_not_nan(const char *function, const T_y &y, const char *name, T_result *result, const Policy &)
Checks if the variable y is nan.
double pi()
Return the value of pi.
Probability, optimization and sampling library.
Definition: agrad.cpp:6
Template specification of functions in std for Stan.
Definition: agrad.hpp:2306
int isnan(const stan::agrad::var &a)
Checks if the given number is NaN.
Definition: agrad.hpp:2361

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