Automatic Differentiation
 
Loading...
Searching...
No Matches
multi_normal_cholesky_lpdf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_PRIM_PROB_MULTI_NORMAL_CHOLESKY_LPDF_HPP
2#define STAN_MATH_PRIM_PROB_MULTI_NORMAL_CHOLESKY_LPDF_HPP
3
18
19namespace stan {
20namespace math {
21
44template <bool propto, typename T_y, typename T_loc, typename T_covar,
45 require_any_not_vector_vt<is_stan_scalar, T_y, T_loc>* = nullptr,
47 T_y, T_loc, T_covar>* = nullptr>
49 const T_y& y, const T_loc& mu, const T_covar& L) {
50 static constexpr const char* function = "multi_normal_cholesky_lpdf";
51 using T_covar_elem = typename scalar_type<T_covar>::type;
53 using T_partials_return = partials_return_t<T_y, T_loc, T_covar>;
54 using matrix_partials_t
55 = Eigen::Matrix<T_partials_return, Eigen::Dynamic, Eigen::Dynamic>;
56 using T_y_ref = ref_type_t<T_y>;
57 using T_mu_ref = ref_type_t<T_loc>;
58 using T_L_ref = ref_type_t<T_covar>;
59
60 check_consistent_sizes_mvt(function, "y", y, "mu", mu);
61 size_t number_of_y = size_mvt(y);
62 size_t number_of_mu = size_mvt(mu);
63 if (number_of_y == 0 || number_of_mu == 0) {
64 return 0;
65 }
66
67 T_y_ref y_ref = y;
68 T_mu_ref mu_ref = mu;
69 T_L_ref L_ref = L;
70 vector_seq_view<T_y_ref> y_vec(y_ref);
71 vector_seq_view<T_mu_ref> mu_vec(mu_ref);
72 const size_t size_vec = max_size_mvt(y, mu);
73
74 const int size_y = y_vec[0].size();
75 const int size_mu = mu_vec[0].size();
76
77 // check size consistency of all random variables y
78 for (size_t i = 1, size_mvt_y = size_mvt(y); i < size_mvt_y; i++) {
79 check_size_match(function,
80 "Size of one of the vectors of "
81 "the random variable",
82 y_vec[i].size(),
83 "Size of the first vector of the "
84 "random variable",
85 size_y);
86 }
87 // check size consistency of all means mu
88 for (size_t i = 1, size_mvt_mu = size_mvt(mu); i < size_mvt_mu; i++) {
89 check_size_match(function,
90 "Size of one of the vectors of "
91 "the location variable",
92 mu_vec[i].size(),
93 "Size of the first vector of the "
94 "location variable",
95 size_mu);
96 }
97
98 check_size_match(function, "Size of random variable", size_y,
99 "size of location parameter", size_mu);
100 check_size_match(function, "Size of random variable", size_y,
101 "rows of covariance parameter", L.rows());
102 check_size_match(function, "Size of random variable", size_y,
103 "columns of covariance parameter", L.cols());
104
105 for (size_t i = 0; i < size_vec; i++) {
106 check_finite(function, "Location parameter", mu_vec[i]);
107 check_not_nan(function, "Random variable", y_vec[i]);
108 }
109 check_cholesky_factor(function, "Cholesky decomposition of a variance matrix",
110 L_ref);
111
112 if (unlikely(size_y == 0)) {
113 return T_return(0);
114 }
115
116 auto ops_partials = make_partials_propagator(y_ref, mu_ref, L_ref);
117
118 T_partials_return logp(0);
120 logp += NEG_LOG_SQRT_TWO_PI * size_y * size_vec;
121 }
122
124 Eigen::Matrix<T_partials_return, Eigen::Dynamic, Eigen::Dynamic>
125 y_val_minus_mu_val(size_y, size_vec);
126
127 for (size_t i = 0; i < size_vec; i++) {
128 decltype(auto) y_val = as_value_column_vector_or_scalar(y_vec[i]);
129 decltype(auto) mu_val = as_value_column_vector_or_scalar(mu_vec[i]);
130 y_val_minus_mu_val.col(i) = y_val - mu_val;
131 }
132
133 matrix_partials_t half;
134 matrix_partials_t scaled_diff;
135
136 // If the covariance is not autodiff, we can avoid computing a matrix
137 // inverse
139 matrix_partials_t L_val = value_of(L_ref);
140
141 half = mdivide_left_tri<Eigen::Lower>(L_val, y_val_minus_mu_val)
142 .transpose();
143
144 scaled_diff = mdivide_right_tri<Eigen::Lower>(half, L_val).transpose();
145
147 logp -= sum(log(L_val.diagonal())) * size_vec;
148 }
149 } else {
150 matrix_partials_t inv_L_val
151 = mdivide_left_tri<Eigen::Lower>(value_of(L_ref));
152
153 half = (inv_L_val.template triangularView<Eigen::Lower>()
154 * y_val_minus_mu_val)
155 .transpose();
156
157 scaled_diff = (half * inv_L_val.template triangularView<Eigen::Lower>())
158 .transpose();
159
160 logp += sum(log(inv_L_val.diagonal())) * size_vec;
161 partials<2>(ops_partials) -= size_vec * inv_L_val.transpose();
162
163 for (size_t i = 0; i < size_vec; i++) {
164 partials_vec<2>(ops_partials)[i] += scaled_diff.col(i) * half.row(i);
165 }
166 }
167
168 logp -= 0.5 * sum(columns_dot_self(half));
169
170 for (size_t i = 0; i < size_vec; i++) {
172 partials_vec<0>(ops_partials)[i] -= scaled_diff.col(i);
173 }
175 partials_vec<1>(ops_partials)[i] += scaled_diff.col(i);
176 }
177 }
178 }
179
180 return ops_partials.build(logp);
181}
182
203template <bool propto, typename T_y, typename T_loc, typename T_covar,
206 T_y, T_loc, T_covar>* = nullptr>
208 const T_y& y, const T_loc& mu, const T_covar& L) {
209 static constexpr const char* function = "multi_normal_cholesky_lpdf";
210 using T_covar_elem = typename scalar_type<T_covar>::type;
211 using T_return = return_type_t<T_y, T_loc, T_covar>;
212 using T_partials_return = partials_return_t<T_y, T_loc, T_covar>;
213 using matrix_partials_t
214 = Eigen::Matrix<T_partials_return, Eigen::Dynamic, Eigen::Dynamic>;
215 using vector_partials_t = Eigen::Matrix<T_partials_return, Eigen::Dynamic, 1>;
216 using row_vector_partials_t
217 = Eigen::Matrix<T_partials_return, 1, Eigen::Dynamic>;
218 using T_y_ref = ref_type_t<T_y>;
219 using T_mu_ref = ref_type_t<T_loc>;
220 using T_L_ref = ref_type_t<T_covar>;
221
222 T_y_ref y_ref = y;
223 T_mu_ref mu_ref = mu;
224 T_L_ref L_ref = L;
225 decltype(auto) y_val = as_value_column_vector_or_scalar(y_ref);
226 decltype(auto) mu_val = as_value_column_vector_or_scalar(mu_ref);
227
228 const int size_y = y_ref.size();
229 const int size_mu = mu_ref.size();
230
231 check_size_match(function, "Size of random variable", size_y,
232 "size of location parameter", size_mu);
233 check_size_match(function, "Size of random variable", size_y,
234 "rows of covariance parameter", L.rows());
235 check_size_match(function, "Size of random variable", size_y,
236 "columns of covariance parameter", L.cols());
237
238 check_finite(function, "Location parameter", mu_val);
239 check_not_nan(function, "Random variable", y_val);
240
241 if (unlikely(size_y == 0)) {
242 return T_return(0);
243 }
244
245 auto ops_partials = make_partials_propagator(y_ref, mu_ref, L_ref);
246
247 T_partials_return logp(0);
248 if (include_summand<propto>::value) {
249 logp += NEG_LOG_SQRT_TWO_PI * size_y;
250 }
251
252 if (include_summand<propto, T_y, T_loc, T_covar_elem>::value) {
253 row_vector_partials_t half;
254 vector_partials_t scaled_diff;
255 vector_partials_t y_val_minus_mu_val = y_val - mu_val;
256
257 // If the covariance is not autodiff, we can avoid computing a matrix
258 // inverse
259 if (is_constant<T_covar_elem>::value) {
260 matrix_partials_t L_val = value_of(L_ref);
261
262 half = mdivide_left_tri<Eigen::Lower>(L_val, y_val_minus_mu_val)
263 .transpose();
264
265 scaled_diff = mdivide_right_tri<Eigen::Lower>(half, L_val).transpose();
266
267 if (include_summand<propto>::value) {
268 logp -= sum(log(L_val.diagonal()));
269 }
270 } else {
271 matrix_partials_t inv_L_val
272 = mdivide_left_tri<Eigen::Lower>(value_of(L_ref));
273
274 half = (inv_L_val.template triangularView<Eigen::Lower>()
275 * y_val_minus_mu_val.template cast<T_partials_return>())
276 .transpose();
277
278 scaled_diff = (half * inv_L_val.template triangularView<Eigen::Lower>())
279 .transpose();
280
281 logp += sum(log(inv_L_val.diagonal()));
282 edge<2>(ops_partials).partials_
283 += scaled_diff * half - inv_L_val.transpose();
284 }
285
286 logp -= 0.5 * sum(dot_self(half));
287
288 if (!is_constant_all<T_y>::value) {
289 partials<0>(ops_partials) -= scaled_diff;
290 }
291 if (!is_constant_all<T_loc>::value) {
292 partials<1>(ops_partials) += scaled_diff;
293 }
294 }
295
296 return ops_partials.build(logp);
297}
298
299template <typename T_y, typename T_loc, typename T_covar>
301 const T_y& y, const T_loc& mu, const T_covar& L) {
302 return multi_normal_cholesky_lpdf<false>(y, mu, L);
303}
304
305} // namespace math
306} // namespace stan
307#endif
This class provides a low-cost wrapper for situations where you either need an Eigen Vector or RowVec...
#define unlikely(x)
require_all_not_t< is_nonscalar_prim_or_rev_kernel_expression< std::decay_t< Types > >... > require_all_not_nonscalar_prim_or_rev_kernel_expression_t
Require none of the types satisfy is_nonscalar_prim_or_rev_kernel_expression.
auto transpose(Arg &&a)
Transposes a kernel generator expression.
return_type_t< T_y_cl, T_loc_cl, T_covar_cl > multi_normal_cholesky_lpdf(const T_y_cl &y, const T_loc_cl &mu, const T_covar_cl &L)
The log of the multivariate normal density for the given y, mu, and a Cholesky factor L of the varian...
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
size_t size_mvt(const ScalarT &)
Provides the size of a multivariate argument.
Definition size_mvt.hpp:24
typename return_type< Ts... >::type return_type_t
Convenience type for the return type of the specified template parameters.
require_all_t< container_type_check_base< is_vector, value_type_t, TypeCheck, Check >... > require_all_vector_vt
Require all of the types satisfy is_vector.
void check_consistent_sizes_mvt(const char *)
Trivial no input case, this function is a no-op.
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 , .
auto columns_dot_self(const T &a)
Returns the dot product of each column of a matrix with itself.
fvar< T > sum(const std::vector< fvar< T > > &m)
Return the sum of the entries of the specified standard vector.
Definition sum.hpp:22
size_t max_size_mvt(const T1 &x1, const Ts &... xs)
Calculate the size of the largest multivariate input.
void check_finite(const char *function, const char *name, const T_y &y)
Return true if all values in y are finite.
void check_not_nan(const char *function, const char *name, const T_y &y)
Check if y is not NaN.
void check_cholesky_factor(const char *function, const char *name, const Mat &y)
Throw an exception if the specified matrix is not a valid Cholesky factor.
auto dot_self(const T &a)
Returns squared norm of a vector or matrix.
Definition dot_self.hpp:21
void check_size_match(const char *function, const char *name_i, T_size1 i, const char *name_j, T_size2 j)
Check if the provided sizes match.
auto as_value_column_vector_or_scalar(T &&a)
Extract values from input argument and transform to a column vector.
auto make_partials_propagator(Ops &&... ops)
Construct an partials_propagator.
typename ref_type_if< true, T >::type ref_type_t
Definition ref_type.hpp:55
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
Metaprogramming struct to detect whether a given type is constant in the mathematical sense (not the ...
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...
std::decay_t< T > type