1 #ifndef STAN__PROB__DISTRIBUTIONS__UNIVARIATE__DISCRETE__NEG_BINOMIAL_2_HPP
2 #define STAN__PROB__DISTRIBUTIONS__UNIVARIATE__DISCRETE__NEG_BINOMIAL_2_HPP
4 #include <boost/random/negative_binomial_distribution.hpp>
5 #include <boost/random/variate_generator.hpp>
7 #include <boost/math/special_functions/digamma.hpp>
27 template <
bool propto,
29 typename T_location,
typename T_inv_scale>
30 typename return_type<T_location, T_inv_scale>::type
33 const T_inv_scale& phi) {
35 static const char*
function =
"stan::prob::neg_binomial_log(%1%)";
57 "Inverse scale parameter",
77 operands_and_partials(mu, phi);
84 for (
size_t i = 0, size =
length(mu); i <
size; ++i)
89 for (
size_t i = 0, size =
length(phi); i <
size; ++i)
94 for (
size_t i = 0, size =
length(phi); i <
size; ++i)
95 log_phi[i] =
log(phi__[i]);
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]);
106 for (
size_t i = 0; i < len_np; ++i)
107 n_plus_phi[i] = n_vec[i] + phi__[i];
109 for (
size_t i = 0; i <
size; i++) {
111 logp -=
lgamma(n_vec[i] + 1.0);
115 logp -= (n_plus_phi[i])*log_mu_plus_phi[i];
119 logp +=
lgamma(n_plus_phi[i]);
122 operands_and_partials.
d_x1[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]);
131 return operands_and_partials.
to_var(logp);
134 template <
typename T_n,
135 typename T_location,
typename T_inv_scale>
139 const T_location& mu,
140 const T_inv_scale& phi) {
141 return neg_binomial_2_log<false>(n, mu, phi);
146 template <
bool propto,
148 typename T_log_location,
typename T_inv_scale>
151 const T_log_location& eta,
152 const T_inv_scale& phi) {
154 static const char*
function =
"stan::prob::neg_binomial_log(%1%)";
171 check_finite(
function, eta,
"Log location parameter", &logp);
176 "Log location parameter",
177 "Inverse scale parameter",
197 operands_and_partials(eta,phi);
204 for (
size_t i = 0, size =
length(eta); i <
size; ++i)
209 for (
size_t i = 0, size =
length(phi); i <
size; ++i)
215 for (
size_t i = 0, size =
length(phi); i <
size; ++i)
216 log_phi[i] =
log(phi__[i]);
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]);
227 for (
size_t i = 0; i < len_np; ++i)
228 n_plus_phi[i] = n_vec[i] + phi__[i];
230 for (
size_t i = 0; i <
size; i++) {
232 logp -=
lgamma(n_vec[i] + 1.0);
236 logp -= (n_plus_phi[i])*logsumexp_eta_logphi[i];
238 logp += n_vec[i]*eta__[i];
240 logp +=
lgamma(n_plus_phi[i]);
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]);
251 return operands_and_partials.
to_var(logp);
254 template <
typename T_n,
255 typename T_log_location,
typename T_inv_scale>
259 const T_log_location& eta,
260 const T_inv_scale& phi) {
261 return neg_binomial_2_log_log<false>(n,eta,phi);
269 using boost::variate_generator;
270 using boost::random::negative_binomial_distribution;
272 static const char*
function =
"stan::prob::neg_binomial_2_rng(%1%)";
289 using boost::variate_generator;
290 using boost::random::negative_binomial_distribution;
292 static const char*
function =
"stan::prob::neg_binomial_2_log_rng(%1%)";
297 check_finite(
function, eta,
"Log-location parameter", (
double*)0);
fvar< T > log_sum_exp(const fvar< T > &x1, const fvar< T > &x2)
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.
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...
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.
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.
int neg_binomial_2_rng(const double mu, const double phi, RNG &rng)
fvar< T > lgamma(const fvar< T > &x)
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
Metaprogram to determine if a type has a base scalar type that can be assigned to type double...
double value_of(const T x)
Return the value of the specified scalar argument converted to a double value.
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...
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)
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)
fvar< T > multiply_log(const fvar< T > &x1, const fvar< T > &x2)
int size(const std::vector< T > &x)
fvar< T > digamma(const fvar< T > &x)
VectorView< double *, is_vector< T1 >::value, is_constant_struct< T1 >::value > d_x1
fvar< T > log(const fvar< T > &x)
VectorView is a template metaprogram that takes its argument and allows it to be used like a vector...
fvar< T > exp(const fvar< T > &x)
double gamma_rng(const double alpha, const double beta, RNG &rng)