Automatic Differentiation
 
Loading...
Searching...
No Matches
cholesky_decompose.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_MAT_FUN_CHOLESKY_DECOMPOSE_HPP
2#define STAN_MATH_REV_MAT_FUN_CHOLESKY_DECOMPOSE_HPP
3
15#include <algorithm>
16#include <vector>
17
18namespace stan {
19namespace math {
20
21namespace internal {
22template <typename LMat, typename LAMat>
23inline void initialize_return(LMat& L, const LAMat& L_A, vari*& dummy) {
24 for (Eigen::Index j = 0; j < L_A.rows(); ++j) {
25 for (Eigen::Index i = 0; i < L_A.rows(); ++i) {
26 if (j > i) {
27 L.coeffRef(i, j) = dummy;
28 } else {
29 L.coeffRef(i, j) = new vari(L_A.coeffRef(i, j), false);
30 }
31 }
32 }
33}
34
47template <typename T1, typename T2, typename T3>
48inline auto unblocked_cholesky_lambda(T1& L_A, T2& L, T3& A) {
49 return [L_A, L, A]() mutable {
50 const size_t N = A.rows();
51 // Algorithm is in rowmajor so we make the adjoint copy rowmajor
52 Eigen::Matrix<double, -1, -1, Eigen::RowMajor> adjL(L.rows(), L.cols());
53 Eigen::MatrixXd adjA = Eigen::MatrixXd::Zero(L.rows(), L.cols());
54 adjL.template triangularView<Eigen::Lower>() = L.adj();
55 for (int i = N - 1; i >= 0; --i) {
56 for (int j = i; j >= 0; --j) {
57 if (i == j) {
58 adjA.coeffRef(i, j) = 0.5 * adjL.coeff(i, j) / L_A.coeff(i, j);
59 } else {
60 adjA.coeffRef(i, j) = adjL.coeff(i, j) / L_A.coeff(j, j);
61 adjL.coeffRef(j, j)
62 -= adjL.coeff(i, j) * L_A.coeff(i, j) / L_A.coeff(j, j);
63 }
64 for (int k = j - 1; k >= 0; --k) {
65 adjL.coeffRef(i, k) -= adjA.coeff(i, j) * L_A.coeff(j, k);
66 adjL.coeffRef(j, k) -= adjA.coeff(i, j) * L_A.coeff(i, k);
67 }
68 }
69 }
70 A.adj() += adjA;
71 };
72}
73
80template <typename T1, typename T2, typename T3>
81inline auto cholesky_lambda(T1& L_A, T2& L, T3& A) {
82 return [L_A, L, A]() mutable {
83 using Eigen::Lower;
84 using Eigen::StrictlyUpper;
85 using Eigen::Upper;
86 Eigen::MatrixXd L_adj = Eigen::MatrixXd::Zero(L.rows(), L.cols());
87 L_adj.template triangularView<Eigen::Lower>() = L.adj();
88 const int M_ = L_A.rows();
89 int block_size_ = std::max(M_ / 8, 8);
90 block_size_ = std::min(block_size_, 128);
91 for (int k = M_; k > 0; k -= block_size_) {
92 int j = std::max(0, k - block_size_);
93 auto R = L_A.block(j, 0, k - j, j);
94 auto D = L_A.block(j, j, k - j, k - j).eval();
95 auto B = L_A.block(k, 0, M_ - k, j);
96 auto C = L_A.block(k, j, M_ - k, k - j);
97 auto R_adj = L_adj.block(j, 0, k - j, j);
98 auto D_adj = L_adj.block(j, j, k - j, k - j);
99 auto B_adj = L_adj.block(k, 0, M_ - k, j);
100 auto C_adj = L_adj.block(k, j, M_ - k, k - j);
101 D.transposeInPlace();
102 if (C_adj.size() > 0) {
103 C_adj = D.template triangularView<Upper>()
104 .solve(C_adj.transpose())
105 .transpose();
106 B_adj.noalias() -= C_adj * R;
107 D_adj.noalias() -= C_adj.transpose() * C;
108 }
109 D_adj = (D * D_adj.template triangularView<Lower>()).eval();
110 D_adj.template triangularView<StrictlyUpper>()
111 = D_adj.adjoint().template triangularView<StrictlyUpper>();
112 D.template triangularView<Upper>().solveInPlace(D_adj);
113 D.template triangularView<Upper>().solveInPlace(D_adj.transpose());
114 R_adj.noalias() -= C_adj.transpose() * B;
115 R_adj.noalias() -= D_adj.template selfadjointView<Lower>() * R;
116 D_adj.diagonal() *= 0.5;
117 }
118 A.adj().template triangularView<Eigen::Lower>() += L_adj;
119 };
120}
121} // namespace internal
122
134template <typename EigMat, require_eigen_vt<is_var, EigMat>* = nullptr>
135inline auto cholesky_decompose(const EigMat& A) {
136 check_square("cholesky_decompose", "A", A);
137 arena_t<EigMat> arena_A = A;
138 arena_t<Eigen::Matrix<double, -1, -1>> L_A(arena_A.val());
139
140 check_symmetric("cholesky_decompose", "A", A);
141 Eigen::LLT<Eigen::Ref<Eigen::MatrixXd>, Eigen::Lower> L_factor(L_A);
142 check_pos_definite("cholesky_decompose", "m", L_factor);
143
144 L_A.template triangularView<Eigen::StrictlyUpper>().setZero();
145 // looping gradient calcs faster for small matrices compared to
146 // cholesky_block
147 vari* dummy = new vari(0.0, false);
148 arena_t<EigMat> L(L_A.rows(), L_A.cols());
149 if (L_A.rows() <= 35) {
150 internal::initialize_return(L, L_A, dummy);
152 } else {
153 internal::initialize_return(L, L_A, dummy);
155 }
156 return plain_type_t<EigMat>(L);
157}
158
170template <typename T, require_var_matrix_t<T>* = nullptr>
171inline auto cholesky_decompose(const T& A) {
172 check_symmetric("cholesky_decompose", "A", A.val());
174 if (A.rows() <= 35) {
176 } else {
178 }
179 return L;
180}
181
182} // namespace math
183} // namespace stan
184#endif
void check_symmetric(const char *function, const char *name, const matrix_cl< T > &y)
Check if the matrix_cl is symmetric.
void initialize_return(LMat &L, const LAMat &L_A, vari *&dummy)
auto unblocked_cholesky_lambda(T1 &L_A, T2 &L, T3 &A)
Reverse mode differentiation algorithm reference:
auto cholesky_lambda(T1 &L_A, T2 &L, T3 &A)
Reverse mode differentiation algorithm reference:
void check_square(const char *function, const char *name, const T_y &y)
Check if the specified matrix is square.
matrix_cl< double > cholesky_decompose(const matrix_cl< double > &A)
Returns the lower-triangular Cholesky factor (i.e., matrix square root) of the specified square,...
void reverse_pass_callback(F &&functor)
Puts a callback on the autodiff stack to be called in reverse pass.
T eval(T &&arg)
Inputs which have a plain_type equal to the own time are forwarded unmodified (for Eigen expressions ...
Definition eval.hpp:20
void check_pos_definite(const char *function, const char *name, const EigMat &y)
Check if the specified square, symmetric matrix is positive definite.
vari_value< double > vari
Definition vari.hpp:197
typename plain_type< T >::type plain_type_t
typename internal::arena_type_impl< std::decay_t< T > >::type arena_t
Determines a type that can be used in place of T that does any dynamic allocations on the AD stack.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...