Stan Math Library  2.20.0
reverse mode automatic differentiation
bernoulli_logit_glm_lpmf.hpp
Go to the documentation of this file.
1 #ifndef STAN_MATH_PRIM_MAT_PROB_BERNOULLI_LOGIT_GLM_LPMF_HPP
2 #define STAN_MATH_PRIM_MAT_PROB_BERNOULLI_LOGIT_GLM_LPMF_HPP
3 
12 #include <cmath>
13 
14 namespace stan {
15 namespace math {
16 
44 template <bool propto, typename T_y, typename T_x, typename T_alpha,
45  typename T_beta>
47  const T_y &y, const T_x &x, const T_alpha &alpha, const T_beta &beta) {
48  static const char *function = "bernoulli_logit_glm_lpmf";
50  T_partials_return;
51  typedef typename std::conditional<
53  Eigen::Matrix<typename partials_return_type<T_y>::type, -1, 1>,
54  typename partials_return_type<T_y>::type>::type T_y_val;
55 
56  using Eigen::Dynamic;
57  using Eigen::Matrix;
58  using std::exp;
59 
60  const size_t N = x.rows();
61  const size_t M = x.cols();
62 
63  check_bounded(function, "Vector of dependent variables", y, 0, 1);
64  check_consistent_size(function, "Vector of dependent variables", y, N);
65  check_consistent_size(function, "Weight vector", beta, M);
67  check_consistent_sizes(function, "Vector of intercepts", alpha,
68  "Vector of dependent variables", y);
69 
70  if (size_zero(y, x, beta))
71  return 0;
72 
74  return 0;
75 
76  T_partials_return logp(0);
77  const auto &x_val = value_of_rec(x);
78  const auto &y_val = value_of_rec(y);
79  const auto &beta_val = value_of_rec(beta);
80  const auto &alpha_val = value_of_rec(alpha);
81 
82  const auto &y_val_vec = as_column_vector_or_scalar(y_val);
83  const auto &beta_val_vec = as_column_vector_or_scalar(beta_val);
84  const auto &alpha_val_vec = as_column_vector_or_scalar(alpha_val);
85 
86  T_y_val signs = 2 * as_array_or_scalar(y_val_vec) - 1;
87 
88  Eigen::Array<T_partials_return, Dynamic, 1> ytheta = x_val * beta_val_vec;
89  ytheta = as_array_or_scalar(signs)
90  * (ytheta + as_array_or_scalar(alpha_val_vec));
91 
92  // Compute the log-density and handle extreme values gracefully
93  // using Taylor approximations.
94  // And compute the derivatives wrt theta.
95  static const double cutoff = 20.0;
96  Eigen::Array<T_partials_return, Dynamic, 1> exp_m_ytheta = exp(-ytheta);
97  logp += sum(
98  (ytheta > cutoff)
99  .select(-exp_m_ytheta,
100  (ytheta < -cutoff).select(ytheta, -log1p(exp_m_ytheta))));
101  if (!std::isfinite(logp)) {
102  check_finite(function, "Weight vector", beta);
103  check_finite(function, "Intercept", alpha);
104  check_finite(function, "Matrix of independent variables", ytheta);
105  }
106 
107  // Compute the necessary derivatives.
108  operands_and_partials<T_x, T_alpha, T_beta> ops_partials(x, alpha, beta);
110  Matrix<T_partials_return, Dynamic, 1> theta_derivative
111  = (ytheta > cutoff)
112  .select(-exp_m_ytheta,
113  (ytheta < -cutoff)
114  .select(as_array_or_scalar(signs),
115  as_array_or_scalar(signs) * exp_m_ytheta
116  / (exp_m_ytheta + 1)));
118  ops_partials.edge3_.partials_ = x_val.transpose() * theta_derivative;
119  }
121  ops_partials.edge1_.partials_
122  = (beta_val_vec * theta_derivative.transpose()).transpose();
123  }
126  ops_partials.edge2_.partials_ = theta_derivative;
127  else
128  ops_partials.edge2_.partials_[0] = sum(theta_derivative);
129  }
130  }
131  return ops_partials.build(logp);
132 }
133 
134 template <typename T_y, typename T_x, typename T_alpha, typename T_beta>
136 bernoulli_logit_glm_lpmf(const T_y &y, const T_x &x, const T_alpha &alpha,
137  const T_beta &beta) {
138  return bernoulli_logit_glm_lpmf<false>(y, x, alpha, beta);
139 }
140 } // namespace math
141 } // namespace stan
142 #endif
fvar< T > sum(const std::vector< fvar< T > > &m)
Return the sum of the entries of the specified standard vector.
Definition: sum.hpp:20
const Eigen::Matrix< T, Eigen::Dynamic, 1 > & as_column_vector_or_scalar(const Eigen::Matrix< T, Eigen::Dynamic, 1 > &a)
Converts input argument to a column vector or a scalar.
void check_finite(const char *function, const char *name, const T_y &y)
Check if y is finite.
bool isfinite(const stan::math::var &v)
Checks if the given number has finite value.
boost::math::tools::promote_args< double, typename partials_type< typename scalar_type< T >::type >::type, typename partials_return_type< T_pack... >::type >::type type
void check_bounded(const char *function, const char *name, const T_y &y, const T_low &low, const T_high &high)
Check if the value is between the low and high values, inclusively.
double value_of_rec(const fvar< T > &v)
Return the value of the specified variable.
Extends std::true_type when instantiated with zero or more template parameters, all of which extend t...
Definition: conjunction.hpp:14
This template builds partial derivatives with respect to a set of operands.
bool size_zero(T &x)
Returns 1 if input is of length 0, returns 0 otherwise.
Definition: size_zero.hpp:18
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
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
void check_consistent_size(const char *function, const char *name, const T &x, size_t expected_size)
Check if the dimension of x is consistent, which is defined to be expected_size if x is a vector or 1...
boost::math::tools::promote_args< double, typename scalar_type< T >::type, typename return_type< Types_pack... >::type >::type type
Definition: return_type.hpp:36
fvar< T > exp(const fvar< T > &x)
Definition: exp.hpp:11
return_type< T_x, T_alpha, T_beta >::type bernoulli_logit_glm_lpmf(const T_y &y, const T_x &x, const T_alpha &alpha, const T_beta &beta)
Returns the log PMF of the Generalized Linear Model (GLM) with Bernoulli distribution and logit link ...
T_return_type build(double value)
Build the node to be stored on the autodiff graph.
matrix_cl transpose(const matrix_cl &src)
Takes the transpose of the matrix on the OpenCL device.
Definition: transpose.hpp:20
internal::ops_partials_edge< double, Op2 > edge2_
Eigen::ArrayWrapper< const Eigen::Matrix< T, R, C > > as_array_or_scalar(const Eigen::Matrix< T, R, C > &v)
Converts a matrix type to an array.
fvar< T > log1p(const fvar< T > &x)
Definition: log1p.hpp:12
void check_consistent_sizes(const char *function, const char *name1, const T1 &x1, const char *name2, const T2 &x2)
Check if the dimension of x1 is consistent with x2.
internal::ops_partials_edge< double, Op3 > edge3_
internal::ops_partials_edge< double, Op1 > edge1_

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