Automatic Differentiation
 
Loading...
Searching...
No Matches
eigendecompose_sym.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUN_EIGENDECOMPOSE_HPP
2#define STAN_MATH_REV_FUN_EIGENDECOMPOSE_HPP
3
13
14namespace stan {
15namespace math {
16
27template <typename T, require_rev_matrix_t<T>* = nullptr>
28inline auto eigendecompose_sym(const T& m) {
29 using eigval_return_t = return_var_matrix_t<Eigen::VectorXd, T>;
30 using eigvec_return_t = return_var_matrix_t<T>;
31
32 if (unlikely(m.size() == 0)) {
33 return std::make_tuple(eigvec_return_t(Eigen::MatrixXd(0, 0)),
34 eigval_return_t(Eigen::VectorXd(0)));
35 }
36 check_symmetric("eigendecompose_sym", "m", m);
37
38 auto arena_m = to_arena(m);
39 Eigen::SelfAdjointEigenSolver<Eigen::MatrixXd> solver(arena_m.val());
40 arena_t<eigval_return_t> eigenvals = solver.eigenvalues();
41 arena_t<eigvec_return_t> eigenvecs = solver.eigenvectors();
42
43 reverse_pass_callback([eigenvals, arena_m, eigenvecs]() mutable {
44 // eigenvalue reverse calculation
45 auto value_adj = eigenvecs.val_op() * eigenvals.adj().asDiagonal()
46 * eigenvecs.val_op().transpose();
47 // eigenvector reverse calculation
48 const auto p = arena_m.val().cols();
49 Eigen::MatrixXd f
50 = (1
51 / (eigenvals.val_op().rowwise().replicate(p).transpose()
52 - eigenvals.val_op().rowwise().replicate(p))
53 .array());
54 f.diagonal().setZero();
55 auto vector_adj
56 = eigenvecs.val_op()
57 * f.cwiseProduct(eigenvecs.val_op().transpose() * eigenvecs.adj_op())
58 * eigenvecs.val_op().transpose();
59
60 arena_m.adj() += value_adj + vector_adj;
61 });
62
63 return std::make_tuple(std::move(eigvec_return_t(eigenvecs)),
64 std::move(eigval_return_t(eigenvals)));
65}
66
67} // namespace math
68} // namespace stan
69#endif
#define unlikely(x)
void check_symmetric(const char *function, const char *name, const matrix_cl< T > &y)
Check if the matrix_cl is symmetric.
void reverse_pass_callback(F &&functor)
Puts a callback on the autodiff stack to be called in reverse pass.
arena_t< T > to_arena(const T &a)
Converts given argument into a type that either has any dynamic allocation on AD stack or schedules i...
Definition to_arena.hpp:25
std::tuple< Eigen::Matrix< value_type_t< EigMat >, -1, -1 >, Eigen::Matrix< value_type_t< EigMat >, -1, 1 > > eigendecompose_sym(const EigMat &m)
Return the eigendecomposition of the specified symmetric matrix.
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.
std::conditional_t< is_any_var_matrix< ReturnType, Types... >::value, stan::math::var_value< stan::math::promote_scalar_t< double, plain_type_t< ReturnType > > >, stan::math::promote_scalar_t< stan::math::var_value< double >, plain_type_t< ReturnType > > > return_var_matrix_t
Given an Eigen type and several inputs, determine if a matrix should be var<Matrix> or Matrix<var>.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...