1 #ifndef __STAN__PROB__DISTRIBUTIONS__UNIVARIATE__CONTINUOUS__INV_GAMMA_HPP__
2 #define __STAN__PROB__DISTRIBUTIONS__UNIVARIATE__CONTINUOUS__INV_GAMMA_HPP__
31 template <
bool propto,
32 typename T_y,
typename T_shape,
typename T_scale,
37 static const char*
function =
"stan::prob::inv_gamma_log(%1%)";
43 using boost::math::tools::promote_args;
56 if (!
check_not_nan(
function, y,
"Random variable", &logp, Policy()))
72 "Random variable",
"Shape parameter",
"Scale parameter",
85 for (
size_t n = 0; n <
length(y); n++) {
86 const double y_dbl =
value_of(y_vec[n]);
96 using boost::math::digamma;
102 for(
size_t n = 0; n <
length(y); n++) {
107 inv_y[n] = 1.0 /
value_of(y_vec[n]);
111 lgamma_alpha(
length(alpha));
113 digamma_alpha(
length(alpha));
114 for (
size_t n = 0; n <
length(alpha); n++) {
118 digamma_alpha[n] = digamma(
value_of(alpha_vec[n]));
124 for (
size_t n = 0; n <
length(beta); n++)
127 for (
size_t n = 0; n < N; n++) {
129 const double alpha_dbl =
value_of(alpha_vec[n]);
130 const double beta_dbl =
value_of(beta_vec[n]);
133 logp -= lgamma_alpha[n];
135 logp += alpha_dbl * log_beta[n];
137 logp -= (alpha_dbl+1.0) * log_y[n];
139 logp -= beta_dbl * inv_y[n];
143 operands_and_partials.
d_x1[n] += -(alpha_dbl+1) * inv_y[n] + beta_dbl * inv_y[n] * inv_y[n];
145 operands_and_partials.
d_x2[n] += -digamma_alpha[n] + log_beta[n] - log_y[n];
147 operands_and_partials.
d_x3[n] += alpha_dbl / beta_dbl - inv_y[n];
149 return operands_and_partials.
to_var(logp);
152 template <
bool propto,
153 typename T_y,
typename T_shape,
typename T_scale>
160 template <
typename T_y,
typename T_shape,
typename T_scale,
166 return inv_gamma_log<false>(y,alpha,beta,Policy());
169 template <
typename T_y,
typename T_shape,
typename T_scale>
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 log(const var &a)
Return the natural log of the specified variable (cmath).
boost::math::tools::promote_args< T_a, T_b >::type multiply_log(T_a a, T_b b)
double value_of(T x)
Return the value of the specified scalar argument converted to a double value.
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.
bool check_consistent_sizes(const char *function, const T1 &x1, const T2 &x2, const char *name1, const char *name2, T_result *result, const Policy &)
bool check_positive(const char *function, const T_y &y, const char *name, T_result *result, const Policy &)
bool check_finite(const char *function, const T_y &y, const char *name, T_result *result, const Policy &)
Checks if the variable y is finite.
boost::math::policies::policy default_policy
Default error-handling policy from Boost.
return_type< T_y, T_shape, T_scale >::type inv_gamma_log(const T_y &y, const T_shape &alpha, const T_scale &beta, const Policy &)
The log of an inverse gamma density for y with the specified shape and scale parameters.
Probability, optimization and sampling library.
size_t length(const T &x)
size_t max_size(const T1 &x1, const T2 &x2)
A variable implementation that stores operands and derivatives with respect to the variable.
VectorView< double *, is_vector< T1 >::value > d_x1
VectorView< double *, is_vector< T2 >::value > d_x2
T_return_type to_var(double logp)
VectorView< double *, is_vector< T3 >::value > d_x3
Metaprogram to determine if a type has a base scalar type that can be assigned to type double.
Metaprogramming struct to detect whether a given type is constant in the mathematical sense (not the ...
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
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