Automatic Differentiation
 
Loading...
Searching...
No Matches
normal_id_glm_lpdf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_PRIM_PROB_NORMAL_ID_GLM_LPDF_HPP
2#define STAN_MATH_PRIM_PROB_NORMAL_ID_GLM_LPDF_HPP
3
17#include <cmath>
18
19namespace stan {
20namespace math {
21
54template <bool propto, typename T_y, typename T_x, typename T_alpha,
55 typename T_beta, typename T_scale, require_matrix_t<T_x>* = nullptr>
57 const T_y& y, const T_x& x, const T_alpha& alpha, const T_beta& beta,
58 const T_scale& sigma) {
59 using Eigen::Array;
60 using Eigen::Dynamic;
61 using Eigen::Matrix;
62 using Eigen::VectorXd;
63 using std::isfinite;
64 constexpr int T_x_rows = T_x::RowsAtCompileTime;
65 using T_partials_return
67 using T_scale_val = typename std::conditional_t<
69 Eigen::Array<partials_return_t<T_scale>, -1, 1>,
71 using T_y_scaled_tmp =
72 typename std::conditional_t<T_x_rows == 1, T_partials_return,
73 Array<T_partials_return, Dynamic, 1>>;
74 using T_y_ref = ref_type_if_not_constant_t<T_y>;
75 using T_x_ref = ref_type_if_not_constant_t<T_x>;
76 using T_alpha_ref = ref_type_if_not_constant_t<T_alpha>;
77 using T_beta_ref = ref_type_if_not_constant_t<T_beta>;
78 using T_sigma_ref = ref_type_if_not_constant_t<T_scale>;
79
80 const size_t N_instances = T_x_rows == 1 ? stan::math::size(y) : x.rows();
81 const size_t N_attributes = x.cols();
82
83 static constexpr const char* function = "normal_id_glm_lpdf";
84 check_consistent_size(function, "Vector of dependent variables", y,
85 N_instances);
86 check_consistent_size(function, "Weight vector", beta, N_attributes);
87 check_consistent_size(function, "Vector of scale parameters", sigma,
88 N_instances);
89 check_consistent_size(function, "Vector of intercepts", alpha, N_instances);
90 T_sigma_ref sigma_ref = sigma;
91 const auto& sigma_val = value_of(sigma_ref);
92 const auto& sigma_val_vec = to_ref(as_column_vector_or_scalar(sigma_val));
93 check_positive_finite(function, "Scale vector", sigma_val_vec);
94
95 if (size_zero(y, sigma)) {
96 return 0;
97 }
99 return 0;
100 }
101
102 T_y_ref y_ref = y;
103 T_x_ref x_ref = x;
104 T_alpha_ref alpha_ref = alpha;
105 T_beta_ref beta_ref = beta;
106
107 const auto& y_val = value_of(y_ref);
108 const auto& x_val = to_ref_if<!is_constant<T_beta>::value>(value_of(x_ref));
109 const auto& alpha_val = value_of(alpha_ref);
110 const auto& beta_val = value_of(beta_ref);
111
112 const auto& y_val_vec = as_column_vector_or_scalar(y_val);
113 const auto& alpha_val_vec = as_column_vector_or_scalar(alpha_val);
114 const auto& beta_val_vec = to_ref_if<!is_constant<T_x>::value>(
116
117 T_scale_val inv_sigma = 1.0 / as_array_or_scalar(sigma_val_vec);
118
119 // the most efficient way to calculate this depends on template parameters
120 T_partials_return y_scaled_sq_sum;
121
122 Array<T_partials_return, Dynamic, 1> y_scaled(N_instances);
123 if (T_x_rows == 1) {
124 T_y_scaled_tmp y_scaled_tmp
125 = forward_as<T_y_scaled_tmp>((x_val * beta_val_vec).coeff(0, 0));
126 y_scaled = (as_array_or_scalar(y_val_vec) - y_scaled_tmp
127 - as_array_or_scalar(alpha_val_vec))
128 * inv_sigma;
129 } else {
130 y_scaled = x_val * beta_val_vec;
131 y_scaled = (as_array_or_scalar(y_val_vec) - y_scaled
132 - as_array_or_scalar(alpha_val_vec))
133 * inv_sigma;
134 }
135
136 auto ops_partials
137 = make_partials_propagator(y_ref, x_ref, alpha_ref, beta_ref, sigma_ref);
138
140 Matrix<T_partials_return, Dynamic, 1> mu_derivative = inv_sigma * y_scaled;
143 partials<0>(ops_partials) = -mu_derivative;
144 } else {
145 partials<0>(ops_partials)[0] = -mu_derivative.sum();
146 }
147 }
149 if (T_x_rows == 1) {
150 edge<1>(ops_partials).partials_
151 = forward_as<Array<T_partials_return, Dynamic, T_x_rows>>(
152 beta_val_vec * sum(mu_derivative));
153 } else {
154 edge<1>(ops_partials).partials_
155 = (beta_val_vec * mu_derivative.transpose()).transpose();
156 }
157 }
159 if (T_x_rows == 1) {
160 edge<3>(ops_partials).partials_
161 = forward_as<Matrix<T_partials_return, 1, Dynamic>>(
162 mu_derivative.sum() * x_val);
163 } else {
164 partials<3>(ops_partials) = mu_derivative.transpose() * x_val;
165 }
166 }
169 partials<2>(ops_partials) = mu_derivative;
170 } else {
171 partials<2>(ops_partials)[0] = sum(mu_derivative);
172 }
173 }
176 Array<T_partials_return, Dynamic, 1> y_scaled_sq = y_scaled * y_scaled;
177 y_scaled_sq_sum = sum(y_scaled_sq);
178 partials<4>(ops_partials) = (y_scaled_sq - 1) * inv_sigma;
179 } else {
180 y_scaled_sq_sum = sum(y_scaled * y_scaled);
181 partials<4>(ops_partials)[0]
182 = (y_scaled_sq_sum - N_instances)
184 }
185 } else {
186 y_scaled_sq_sum = sum(y_scaled * y_scaled);
187 }
188 } else {
189 y_scaled_sq_sum = sum(y_scaled * y_scaled);
190 }
191
192 if (!isfinite(y_scaled_sq_sum)) {
193 check_finite(function, "Vector of dependent variables", y_val_vec);
194 check_finite(function, "Weight vector", beta_val_vec);
195 check_finite(function, "Intercept", alpha_val_vec);
196 // if all other checks passed, next will only fail if x is not finite
197 check_finite(function, "Matrix of independent variables", y_scaled_sq_sum);
198 }
199
200 // Compute log probability.
201 T_partials_return logp(0.0);
203 logp += NEG_LOG_SQRT_TWO_PI * N_instances;
204 }
207 logp -= sum(log(sigma_val_vec));
208 } else {
209 logp -= N_instances
211 }
212 }
213 logp -= 0.5 * y_scaled_sq_sum;
214
215 return ops_partials.build(logp);
216}
217
218template <typename T_y, typename T_x, typename T_alpha, typename T_beta,
219 typename T_scale>
221 const T_y& y, const T_x& x, const T_alpha& alpha, const T_beta& beta,
222 const T_scale& sigma) {
223 return normal_id_glm_lpdf<false>(y, x, alpha, beta, sigma);
224}
225} // namespace math
226} // namespace stan
227#endif
isfinite_< as_operation_cl_t< T > > isfinite(T &&a)
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_y_cl, T_x_cl, T_alpha_cl, T_beta_cl, T_sigma_cl > normal_id_glm_lpdf(const T_y_cl &y, const T_x_cl &x, const T_alpha_cl &alpha, const T_beta_cl &beta, const T_sigma_cl &sigma)
Returns the log PDF of the Generalized Linear Model (GLM) with Normal distribution and id link functi...
size_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:18
T_actual && forward_as(T_actual &&a)
Assume which type we get.
typename return_type< Ts... >::type return_type_t
Convenience type for the return type of the specified template parameters.
T as_array_or_scalar(T &&v)
Returns specified input value.
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:15
const double NEG_LOG_SQRT_TWO_PI
The value of minus the natural logarithm of the square root of , .
fvar< T > sum(const std::vector< fvar< T > > &m)
Return the sum of the entries of the specified standard vector.
Definition sum.hpp:22
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
Definition to_ref.hpp:17
void check_finite(const char *function, const char *name, const T_y &y)
Return true if all values in y are finite.
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.
typename ref_type_if<!is_constant< T >::value, T >::type ref_type_if_not_constant_t
Definition ref_type.hpp:62
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 ...
Definition fvar.hpp:9
If the input type T is either an eigen matrix with 1 column or 1 row at compile time or a standard ve...
Extends std::true_type when instantiated with zero or more template parameters, all of which extend t...
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...