Automatic Differentiation
 
Loading...
Searching...
No Matches
neg_binomial_2_log_glm_lpmf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_PRIM_PROB_NEG_BINOMIAL_2_LOG_GLM_LPMF_HPP
2#define STAN_MATH_PRIM_PROB_NEG_BINOMIAL_2_LOG_GLM_LPMF_HPP
3
21#include <vector>
22#include <cmath>
23
24namespace stan {
25namespace math {
26
64template <bool propto, typename T_y, typename T_x, typename T_alpha,
65 typename T_beta, typename T_precision,
66 require_matrix_t<T_x>* = nullptr>
67inline return_type_t<T_x, T_alpha, T_beta, T_precision>
68neg_binomial_2_log_glm_lpmf(const T_y& y, const T_x& x, const T_alpha& alpha,
69 const T_beta& beta, const T_precision& phi) {
70 using Eigen::Array;
71 using Eigen::Dynamic;
72 using Eigen::exp;
73 using Eigen::log1p;
74 using Eigen::Matrix;
75 constexpr int T_x_rows = T_x::RowsAtCompileTime;
76 using T_partials_return
78 using T_precision_val = typename std::conditional_t<
80 Eigen::Array<partials_return_t<T_precision>, -1, 1>,
82 using T_sum_val = typename std::conditional_t<
84 Eigen::Array<partials_return_t<T_y, T_precision>, -1, 1>,
86 using T_theta_tmp =
87 typename std::conditional_t<T_x_rows == 1, T_partials_return,
88 Array<T_partials_return, Dynamic, 1>>;
89 using T_x_ref = ref_type_if_not_constant_t<T_x>;
90 using T_alpha_ref = ref_type_if_not_constant_t<T_alpha>;
91 using T_beta_ref = ref_type_if_not_constant_t<T_beta>;
93
94 const size_t N_instances = T_x_rows == 1 ? stan::math::size(y) : x.rows();
95 const size_t N_attributes = x.cols();
96
97 static constexpr const char* function = "neg_binomial_2_log_glm_lpmf";
98 check_consistent_size(function, "Vector of dependent variables", y,
99 N_instances);
100 check_consistent_size(function, "Weight vector", beta, N_attributes);
101 check_consistent_size(function, "Vector of precision parameters", phi,
102 N_instances);
103 check_consistent_size(function, "Vector of intercepts", alpha, N_instances);
104 T_alpha_ref alpha_ref = alpha;
105 T_beta_ref beta_ref = beta;
106 const auto& beta_val = value_of(beta_ref);
107 const auto& alpha_val = value_of(alpha_ref);
108 const auto& beta_val_vec = to_ref(as_column_vector_or_scalar(beta_val));
109 const auto& alpha_val_vec = to_ref(as_column_vector_or_scalar(alpha_val));
110 check_finite(function, "Weight vector", beta_val_vec);
111 check_finite(function, "Intercept", alpha_val_vec);
112
113 if (size_zero(y, phi)) {
114 return 0;
115 }
116
117 const auto& y_ref = to_ref(y);
118 T_phi_ref phi_ref = phi;
119
120 const auto& y_val = value_of(y_ref);
121 const auto& phi_val = value_of(phi_ref);
122
123 const auto& y_val_vec = to_ref(as_column_vector_or_scalar(y_val));
124 const auto& phi_val_vec = to_ref(as_column_vector_or_scalar(phi_val));
125 check_nonnegative(function, "Failures variables", y_val_vec);
126 check_positive_finite(function, "Precision parameter", phi_val_vec);
127
128 if constexpr (!include_summand<propto, T_x, T_alpha, T_beta,
129 T_precision>::value) {
130 return 0;
131 }
132
133 T_x_ref x_ref = x;
134
135 const auto& x_val = to_ref_if<is_autodiff_v<T_beta>>(value_of(x_ref));
136
137 const auto& y_arr = as_array_or_scalar(y_val_vec);
138 const auto& phi_arr = as_array_or_scalar(phi_val_vec);
139
140 Array<T_partials_return, Dynamic, 1> theta(N_instances);
141 if constexpr (T_x_rows == 1) {
142 T_theta_tmp theta_tmp = (x_val * beta_val_vec)(0, 0);
143 theta = theta_tmp + as_array_or_scalar(alpha_val_vec);
144 } else {
145 theta = (x_val * beta_val_vec).array();
146 theta += as_array_or_scalar(alpha_val_vec);
147 }
148 check_finite(function, "Matrix of independent variables", theta);
149 T_precision_val log_phi = log(phi_arr);
150 Array<T_partials_return, Dynamic, 1> logsumexp_theta_logphi
151 = (theta > log_phi)
152 .select(theta + log1p_exp(log_phi - theta),
153 log_phi + log1p_exp(theta - log_phi));
154
155 T_sum_val y_plus_phi = y_arr + phi_arr;
156
157 // Compute the log-density.
158 T_partials_return logp(0);
159 if constexpr (include_summand<propto>::value) {
160 if constexpr (is_vector<T_y>::value) {
161 logp -= sum(lgamma(y_arr + 1.0));
162 } else {
163 logp -= sum(lgamma(y_arr + 1.0)) * N_instances;
164 }
165 }
167 if constexpr (is_vector<T_precision>::value) {
168 scalar_seq_view<decltype(phi_val_vec)> phi_vec(phi_val_vec);
169 for (size_t n = 0; n < N_instances; ++n) {
170 logp += multiply_log(phi_vec[n], phi_vec[n]) - lgamma(phi_vec[n]);
171 }
172 } else {
173 logp += N_instances * (multiply_log(phi_val, phi_val) - lgamma(phi_val));
174 }
175 }
176 logp -= sum(y_plus_phi * logsumexp_theta_logphi);
177
179 logp += sum(y_arr * theta);
180 }
183 logp += sum(lgamma(y_plus_phi));
184 } else {
185 logp += sum(lgamma(y_plus_phi)) * N_instances;
186 }
187 }
188
189 // Compute the necessary derivatives.
190 auto ops_partials
191 = make_partials_propagator(x_ref, alpha_ref, beta_ref, phi_ref);
192 if constexpr (is_any_autodiff_v<T_x, T_beta, T_alpha, T_precision>) {
193 Array<T_partials_return, Dynamic, 1> theta_exp = theta.exp();
194 if constexpr (is_any_autodiff_v<T_x, T_beta, T_alpha>) {
195 Matrix<T_partials_return, Dynamic, 1> theta_derivative
196 = y_arr - theta_exp * y_plus_phi / (theta_exp + phi_arr);
197 if constexpr (is_autodiff_v<T_beta>) {
198 if constexpr (T_x_rows == 1) {
199 edge<2>(ops_partials).partials_ = theta_derivative.sum() * x_val;
200 } else {
201 edge<2>(ops_partials).partials_
202 = x_val.transpose() * theta_derivative;
203 }
204 }
205 if constexpr (is_autodiff_v<T_x>) {
206 if constexpr (T_x_rows == 1) {
207 edge<0>(ops_partials).partials_
208 = beta_val_vec * theta_derivative.sum();
209 } else {
210 edge<0>(ops_partials).partials_
211 = (beta_val_vec * theta_derivative.transpose()).transpose();
212 }
213 }
214 if constexpr (is_autodiff_v<T_alpha>) {
215 if constexpr (is_vector<T_alpha>::value) {
216 partials<1>(ops_partials) = std::move(theta_derivative);
217 } else {
218 partials<1>(ops_partials)[0] = sum(theta_derivative);
219 }
220 }
221 }
222 if constexpr (is_autodiff_v<T_precision>) {
223 if constexpr (is_vector<T_precision>::value) {
224 edge<3>(ops_partials).partials_
225 = 1 - y_plus_phi / (theta_exp + phi_arr) + log_phi
226 - logsumexp_theta_logphi + digamma(y_plus_phi) - digamma(phi_arr);
227 } else {
228 partials<3>(ops_partials)[0]
229 = N_instances
230 + sum(-y_plus_phi / (theta_exp + phi_arr) + log_phi
231 - logsumexp_theta_logphi + digamma(y_plus_phi)
232 - digamma(phi_arr));
233 }
234 }
235 }
236 return ops_partials.build(logp);
237}
238
239template <typename T_y, typename T_x, typename T_alpha, typename T_beta,
240 typename T_precision>
242neg_binomial_2_log_glm_lpmf(const T_y& y, const T_x& x, const T_alpha& alpha,
243 const T_beta& beta, const T_precision& phi) {
244 return neg_binomial_2_log_glm_lpmf<false>(y, x, alpha, beta, phi);
245}
246} // namespace math
247} // namespace stan
248#endif
scalar_seq_view provides a uniform sequence-like wrapper around either a scalar or a sequence of scal...
select_< as_operation_cl_t< T_condition >, as_operation_cl_t< T_then >, as_operation_cl_t< T_else > > select(T_condition &&condition, T_then &&then, T_else &&els)
Selection operation on kernel generator expressions.
Definition select.hpp:148
auto as_column_vector_or_scalar(T &&a)
as_column_vector_or_scalar of a kernel generator expression.
auto transpose(Arg &&a)
Transposes a kernel generator expression.
return_type_t< T_x_cl, T_alpha_cl, T_beta_cl, T_phi_cl > neg_binomial_2_log_glm_lpmf(const T_y_cl &y, const T_x_cl &x, const T_alpha_cl &alpha, const T_beta_cl &beta, const T_phi_cl &phi)
Returns the log PMF of the Generalized Linear Model (GLM) with Negative-Binomial-2 distribution and l...
typename return_type< Ts... >::type return_type_t
Convenience type for the return type of the specified template parameters.
int64_t size(const T &m)
Returns the size (number of the elements) of a matrix_cl or var_value<matrix_cl<T>>.
Definition size.hpp:19
fvar< T > multiply_log(const fvar< T > &x1, const fvar< T > &x2)
T as_array_or_scalar(T &&v)
Returns specified input value.
void check_nonnegative(const char *function, const char *name, const T_y &y)
Check if y is non-negative.
bool size_zero(const T &x)
Returns 1 if input is of length 0, returns 0 otherwise.
Definition size_zero.hpp:19
void check_consistent_size(const char *function, const char *name, const T &x, size_t expected_size)
Check if x is consistent with size expected_size.
T value_of(const fvar< T > &v)
Return the value of the specified variable.
Definition value_of.hpp:18
fvar< T > log(const fvar< T > &x)
Definition log.hpp:18
fvar< T > log1p_exp(const fvar< T > &x)
Definition log1p_exp.hpp:14
void check_finite(const char *function, const char *name, const T_y &y)
Return true if all values in y are finite.
fvar< T > lgamma(const fvar< T > &x)
Return the natural logarithm of the gamma function applied to the specified argument.
Definition lgamma.hpp:21
auto sum(const std::vector< T > &m)
Return the sum of the entries of the specified standard vector.
Definition sum.hpp:23
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
Definition to_ref.hpp:18
fvar< T > beta(const fvar< T > &x1, const fvar< T > &x2)
Return fvar with the beta function applied to the specified arguments and its gradient.
Definition beta.hpp:51
auto make_partials_propagator(Ops &&... ops)
Construct an partials_propagator.
void check_positive_finite(const char *function, const char *name, const T_y &y)
Check if y is positive and finite.
fvar< T > digamma(const fvar< T > &x)
Return the derivative of the log gamma function at the specified argument.
Definition digamma.hpp:23
typename ref_type_if< is_autodiff_v< T >, T >::type ref_type_if_not_constant_t
Definition ref_type.hpp:63
typename partials_return_type< Args... >::type partials_return_t
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
If the input type T is either an eigen matrix with 1 column or 1 row at compile time or a standard ve...
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...