Automatic Differentiation
 
Loading...
Searching...
No Matches
multi_normal_prec_lpdf.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_PRIM_PROB_MULTI_NORMAL_PREC_LPDF_HPP
2#define STAN_MATH_PRIM_PROB_MULTI_NORMAL_PREC_LPDF_HPP
3
15
16namespace stan {
17namespace math {
18
19template <bool propto, typename T_y, typename T_loc, typename T_covar>
21 const T_y& y, const T_loc& mu, const T_covar& Sigma) {
22 using T_covar_elem = typename scalar_type<T_covar>::type;
24 using Eigen::Matrix;
25 using std::vector;
26 static constexpr const char* function = "multi_normal_prec_lpdf";
27 check_positive(function, "Precision matrix rows", Sigma.rows());
28
29 size_t number_of_y = size_mvt(y);
30 size_t number_of_mu = size_mvt(mu);
31 if (number_of_y == 0 || number_of_mu == 0) {
32 return 0;
33 }
34 check_consistent_sizes_mvt(function, "y", y, "mu", mu);
35
36 lp_type lp(0);
37 vector_seq_view<T_y> y_vec(y);
38 vector_seq_view<T_loc> mu_vec(mu);
39 size_t size_vec = max_size_mvt(y, mu);
40
41 int size_y = y_vec[0].size();
42 int size_mu = mu_vec[0].size();
43 if (size_vec > 1) {
44 for (size_t i = 1, size_mvt_y = size_mvt(y); i < size_mvt_y; i++) {
45 check_size_match(function,
46 "Size of one of the vectors "
47 "of the random variable",
48 y_vec[i].size(),
49 "Size of the first vector of "
50 "the random variable",
51 size_y);
52 }
53 for (size_t i = 1, size_mvt_mu = size_mvt(mu); i < size_mvt_mu; i++) {
54 check_size_match(function,
55 "Size of one of the vectors "
56 "of the location variable",
57 mu_vec[i].size(),
58 "Size of the first vector of "
59 "the location variable",
60 size_mu);
61 }
62 }
63
64 check_size_match(function, "Size of random variable", size_y,
65 "size of location parameter", size_mu);
66 check_size_match(function, "Size of random variable", size_y,
67 "rows of covariance parameter", Sigma.rows());
68 check_size_match(function, "Size of random variable", size_y,
69 "columns of covariance parameter", Sigma.cols());
70
71 for (size_t i = 0; i < size_vec; i++) {
72 check_finite(function, "Location parameter", mu_vec[i]);
73 check_not_nan(function, "Random variable", y_vec[i]);
74 }
75 const auto& Sigma_ref = to_ref(Sigma);
76 check_symmetric(function, "Precision matrix", Sigma_ref);
77
78 auto ldlt_Sigma = make_ldlt_factor(Sigma_ref);
79 check_ldlt_factor(function, "LDLT_Factor of precision parameter", ldlt_Sigma);
80
81 if (size_y == 0) {
82 return lp;
83 }
84
86 lp += 0.5 * log_determinant_ldlt(ldlt_Sigma) * size_vec;
87 }
88
90 lp += NEG_LOG_SQRT_TWO_PI * size_y * size_vec;
91 }
92
94 lp_type sum_lp_vec(0.0);
95 for (size_t i = 0; i < size_vec; i++) {
96 const auto& y_col = as_column_vector_or_scalar(y_vec[i]);
97 const auto& mu_col = as_column_vector_or_scalar(mu_vec[i]);
98 sum_lp_vec += trace_quad_form(Sigma_ref, y_col - mu_col);
99 }
100 lp -= 0.5 * sum_lp_vec;
101 }
102 return lp;
103}
104
105template <typename T_y, typename T_loc, typename T_covar>
107 const T_y& y, const T_loc& mu, const T_covar& Sigma) {
108 return multi_normal_prec_lpdf<false>(y, mu, Sigma);
109}
110
111} // namespace math
112} // namespace stan
113#endif
This class provides a low-cost wrapper for situations where you either need an Eigen Vector or RowVec...
void check_symmetric(const char *function, const char *name, const matrix_cl< T > &y)
Check if the matrix_cl is symmetric.
auto as_column_vector_or_scalar(T &&a)
as_column_vector_or_scalar of a kernel generator expression.
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.
value_type_t< T > log_determinant_ldlt(LDLT_factor< T > &A)
Returns log(abs(det(A))) given a LDLT_factor of A.
auto make_ldlt_factor(const T &A)
Make an LDLT_factor with matrix type plain_type_t<T>
void check_consistent_sizes_mvt(const char *)
Trivial no input case, this function is a no-op.
return_type_t< T_y, T_loc, T_covar > multi_normal_prec_lpdf(const T_y &y, const T_loc &mu, const T_covar &Sigma)
const double NEG_LOG_SQRT_TWO_PI
The value of minus the natural logarithm of the square root of , .
ref_type_t< T && > to_ref(T &&a)
This evaluates expensive Eigen expressions.
Definition to_ref.hpp:17
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_ldlt_factor(const char *function, const char *name, LDLT_factor< T > &A)
Raise domain error if the specified LDLT factor is invalid.
void check_not_nan(const char *function, const char *name, const T_y &y)
Check if y is not NaN.
void check_positive(const char *function, const char *name, const T_y &y)
Check if y is positive.
return_type_t< EigMat1, EigMat2 > trace_quad_form(const EigMat1 &A, const EigMat2 &B)
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.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...
Definition fvar.hpp:9
Template metaprogram to calculate whether a summand needs to be included in a proportional (log) prob...
std::decay_t< T > type