Stan  2.5.0
probability, sampling & optimization
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
neg_binomial_2.hpp
Go to the documentation of this file.
1 #ifndef STAN__PROB__DISTRIBUTIONS__UNIVARIATE__DISCRETE__NEG_BINOMIAL_2_HPP
2 #define STAN__PROB__DISTRIBUTIONS__UNIVARIATE__DISCRETE__NEG_BINOMIAL_2_HPP
3 
4 #include <boost/random/negative_binomial_distribution.hpp>
5 #include <boost/random/variate_generator.hpp>
6 
7 #include <boost/math/special_functions/digamma.hpp>
10 #include <stan/math/constants.hpp>
14 #include <stan/meta/traits.hpp>
15 #include <stan/prob/traits.hpp>
16 #include <stan/prob/constants.hpp>
21 
22 namespace stan {
23 
24  namespace prob {
25 
26  // NegBinomial(n|mu,phi) [mu >= 0; phi > 0; n >= 0]
27  template <bool propto,
28  typename T_n,
29  typename T_location, typename T_inv_scale>
30  typename return_type<T_location, T_inv_scale>::type
31  neg_binomial_2_log(const T_n& n,
32  const T_location& mu,
33  const T_inv_scale& phi) {
34 
35  static const char* function = "stan::prob::neg_binomial_log(%1%)";
36 
42 
43  // check if any vectors are zero length
44  if (!(stan::length(n)
45  && stan::length(mu)
46  && stan::length(phi)))
47  return 0.0;
48 
49  double logp(0.0);
50  check_nonnegative(function, n, "Failures variable", &logp);
51  check_positive_finite(function, mu, "Location parameter", &logp);
52  check_positive_finite(function, phi, "Inverse scale parameter", &logp);
53  check_consistent_sizes(function,
54  n,mu,phi,
55  "Failures variable",
56  "Location parameter",
57  "Inverse scale parameter",
58  &logp);
59 
60  // check if no variables are involved and prop-to
62  return 0.0;
63 
68  using boost::math::lgamma;
69 
70  // set up template expressions wrapping scalars into vector views
71  VectorView<const T_n> n_vec(n);
73  VectorView<const T_inv_scale> phi_vec(phi);
74  size_t size = max_size(n, mu, phi);
75 
77  operands_and_partials(mu, phi);
78 
79  size_t len_ep = max_size(mu, phi);
80  size_t len_np = max_size(n, phi);
81 
83  mu__(length(mu));
84  for (size_t i = 0, size = length(mu); i < size; ++i)
85  mu__[i] = value_of(mu_vec[i]);
86 
88  phi__(length(phi));
89  for (size_t i = 0, size = length(phi); i < size; ++i)
90  phi__[i] = value_of(phi_vec[i]);
91 
93  log_phi(length(phi));
94  for (size_t i = 0, size = length(phi); i < size; ++i)
95  log_phi[i] = log(phi__[i]);
96 
99  log_mu_plus_phi(len_ep);
100  for (size_t i = 0; i < len_ep; ++i)
101  log_mu_plus_phi[i] = log(mu__[i] + phi__[i]);
102 
105  n_plus_phi(len_np);
106  for (size_t i = 0; i < len_np; ++i)
107  n_plus_phi[i] = n_vec[i] + phi__[i];
108 
109  for (size_t i = 0; i < size; i++) {
111  logp -= lgamma(n_vec[i] + 1.0);
113  logp += multiply_log(phi__[i], phi__[i]) - lgamma(phi__[i]);
115  logp -= (n_plus_phi[i])*log_mu_plus_phi[i];
117  logp += multiply_log(n_vec[i], mu__[i]);
119  logp += lgamma(n_plus_phi[i]);
120 
122  operands_and_partials.d_x1[i]
123  += n_vec[i]/mu__[i]
124  - (n_vec[i] + phi__[i])
125  / (mu__[i] + phi__[i]);
127  operands_and_partials.d_x2[i]
128  += 1.0 - n_plus_phi[i]/(mu__[i] + phi__[i])
129  + log_phi[i] - log_mu_plus_phi[i] - digamma(phi__[i]) + digamma(n_plus_phi[i]);
130  }
131  return operands_and_partials.to_var(logp);
132  }
133 
134  template <typename T_n,
135  typename T_location, typename T_inv_scale>
136  inline
138  neg_binomial_2_log(const T_n& n,
139  const T_location& mu,
140  const T_inv_scale& phi) {
141  return neg_binomial_2_log<false>(n, mu, phi);
142  }
143 
144 
145  // NegBinomial(n|eta,phi) [phi > 0; n >= 0]
146  template <bool propto,
147  typename T_n,
148  typename T_log_location, typename T_inv_scale>
151  const T_log_location& eta,
152  const T_inv_scale& phi) {
153 
154  static const char* function = "stan::prob::neg_binomial_log(%1%)";
155 
159  using stan::math::value_of;
162 
163  // check if any vectors are zero length
164  if (!(stan::length(n)
165  && stan::length(eta)
166  && stan::length(phi)))
167  return 0.0;
168 
169  double logp(0.0);
170  check_nonnegative(function, n, "Failures variable", &logp);
171  check_finite(function, eta, "Log location parameter", &logp);
172  check_positive_finite(function, phi, "Inverse scale parameter", &logp);
173  check_consistent_sizes(function,
174  n,eta,phi,
175  "Failures variable",
176  "Log location parameter",
177  "Inverse scale parameter",
178  &logp);
179 
180  // check if no variables are involved and prop-to
182  return 0.0;
183 
187  using boost::math::digamma;
188  using boost::math::lgamma;
189 
190  // set up template expressions wrapping scalars into vector views
191  VectorView<const T_n> n_vec(n);
193  VectorView<const T_inv_scale> phi_vec(phi);
194  size_t size = max_size(n, eta, phi);
195 
197  operands_and_partials(eta,phi);
198 
199  size_t len_ep = max_size(eta, phi);
200  size_t len_np = max_size(n, phi);
201 
203  eta__(length(eta));
204  for (size_t i = 0, size = length(eta); i < size; ++i)
205  eta__[i] = value_of(eta_vec[i]);
206 
208  phi__(length(phi));
209  for (size_t i = 0, size = length(phi); i < size; ++i)
210  phi__[i] = value_of(phi_vec[i]);
211 
212 
214  log_phi(length(phi));
215  for (size_t i = 0, size = length(phi); i < size; ++i)
216  log_phi[i] = log(phi__[i]);
217 
220  logsumexp_eta_logphi(len_ep);
221  for (size_t i = 0; i < len_ep; ++i)
222  logsumexp_eta_logphi[i] = log_sum_exp(eta__[i], log_phi[i]);
223 
226  n_plus_phi(len_np);
227  for (size_t i = 0; i < len_np; ++i)
228  n_plus_phi[i] = n_vec[i] + phi__[i];
229 
230  for (size_t i = 0; i < size; i++) {
232  logp -= lgamma(n_vec[i] + 1.0);
234  logp += multiply_log(phi__[i], phi__[i]) - lgamma(phi__[i]);
236  logp -= (n_plus_phi[i])*logsumexp_eta_logphi[i];
238  logp += n_vec[i]*eta__[i];
240  logp += lgamma(n_plus_phi[i]);
241 
243  operands_and_partials.d_x1[i]
244  += n_vec[i] - n_plus_phi[i]
245  / (phi__[i]/exp(eta__[i]) + 1.0);
247  operands_and_partials.d_x2[i]
248  += 1.0 - n_plus_phi[i]/(exp(eta__[i]) + phi__[i])
249  + log_phi[i] - logsumexp_eta_logphi[i] - digamma(phi__[i]) + digamma(n_plus_phi[i]);
250  }
251  return operands_and_partials.to_var(logp);
252  }
253 
254  template <typename T_n,
255  typename T_log_location, typename T_inv_scale>
256  inline
259  const T_log_location& eta,
260  const T_inv_scale& phi) {
261  return neg_binomial_2_log_log<false>(n,eta,phi);
262  }
263 
264  template <class RNG>
265  inline int
266  neg_binomial_2_rng(const double mu,
267  const double phi,
268  RNG& rng) {
269  using boost::variate_generator;
270  using boost::random::negative_binomial_distribution;
271 
272  static const char* function = "stan::prob::neg_binomial_2_rng(%1%)";
273 
275 
276  check_positive_finite(function, mu, "Location parameter", (double*)0);
277  check_positive_finite(function, phi, "Inverse scale parameter",
278  (double*)0);
279 
281  rng),rng);
282  }
283 
284  template <class RNG>
285  inline int
286  neg_binomial_2_log_rng(const double eta,
287  const double phi,
288  RNG& rng) {
289  using boost::variate_generator;
290  using boost::random::negative_binomial_distribution;
291 
292  static const char* function = "stan::prob::neg_binomial_2_log_rng(%1%)";
293 
296 
297  check_finite(function, eta, "Log-location parameter", (double*)0);
298  check_positive_finite(function, phi, "Inverse scale parameter",
299  (double*)0);
300 
302  rng),rng);
303  }
304  }
305 }
306 #endif
fvar< T > log_sum_exp(const fvar< T > &x1, const fvar< T > &x2)
Definition: log_sum_exp.hpp:15
boost::math::tools::promote_args< T1, T2 >::type log_sum_exp(const T2 &a, const T1 &b)
Calculates the log sum of exponetials without overflow.
Definition: log_sum_exp.hpp:49
T_return_type to_var(double logp)
bool check_positive_finite(const char *function, const T_y &y, const char *name, T_result *result)
int neg_binomial_2_log_rng(const double eta, const double phi, RNG &rng)
boost::math::tools::promote_args< T_a, T_b >::type multiply_log(const T_a a, const T_b b)
Calculated the value of the first argument times log of the second argument while behaving properly w...
size_t length(const T &)
Definition: traits.hpp:159
boost::math::tools::promote_args< T_N, T_n >::type binomial_coefficient_log(const T_N N, const T_n n)
Return the log of the binomial coefficient for the specified arguments.
DoubleVectorView allocates double values to be used as intermediate values.
Definition: traits.hpp:358
bool check_finite(const char *function, const T_y &y, const char *name, T_result *result)
Checks if the variable y is finite.
T value_of(const fvar< T > &v)
Return the value of the specified variable.
Definition: value_of.hpp:16
int neg_binomial_2_rng(const double mu, const double phi, RNG &rng)
fvar< T > lgamma(const fvar< T > &x)
Definition: lgamma.hpp:15
A variable implementation that stores operands and derivatives with respect to the variable...
boost::math::tools::promote_args< typename scalar_type< T1 >::type, typename scalar_type< T2 >::type, typename scalar_type< T3 >::type, typename scalar_type< T4 >::type, typename scalar_type< T5 >::type, typename scalar_type< T6 >::type >::type type
Definition: traits.hpp:406
Metaprogram to determine if a type has a base scalar type that can be assigned to type double...
Definition: traits.hpp:57
double value_of(const T x)
Return the value of the specified scalar argument converted to a double value.
Definition: value_of.hpp:24
return_type< T_log_location, T_inv_scale >::type neg_binomial_2_log_log(const T_n &n, const T_log_location &eta, const T_inv_scale &phi)
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
Definition: traits.hpp:35
return_type< T_location, T_inv_scale >::type neg_binomial_2_log(const T_n &n, const T_location &mu, const T_inv_scale &phi)
int poisson_rng(const double lambda, RNG &rng)
Definition: poisson.hpp:384
VectorView< double *, is_vector< T2 >::value, is_constant_struct< T2 >::value > d_x2
bool check_nonnegative(const char *function, const T_y &y, const char *name, T_result *result)
bool check_consistent_sizes(const char *function, const T1 &x1, const T2 &x2, const char *name1, const char *name2, T_result *result)
size_t max_size(const T1 &x1, const T2 &x2)
Definition: traits.hpp:191
fvar< T > multiply_log(const fvar< T > &x1, const fvar< T > &x2)
int size(const std::vector< T > &x)
Definition: size.hpp:11
fvar< T > digamma(const fvar< T > &x)
Definition: digamma.hpp:16
VectorView< double *, is_vector< T1 >::value, is_constant_struct< T1 >::value > d_x1
fvar< T > log(const fvar< T > &x)
Definition: log.hpp:15
VectorView is a template metaprogram that takes its argument and allows it to be used like a vector...
Definition: traits.hpp:275
fvar< T > exp(const fvar< T > &x)
Definition: exp.hpp:16
double gamma_rng(const double alpha, const double beta, RNG &rng)
Definition: gamma.hpp:491

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