Automatic Differentiation
 
Loading...
Searching...
No Matches
svd_V.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUN_SVD_V_HPP
2#define STAN_MATH_REV_FUN_SVD_V_HPP
3
11
12namespace stan {
13namespace math {
14
25template <typename EigMat, require_rev_matrix_t<EigMat>* = nullptr>
26inline auto svd_V(const EigMat& m) {
28 if (unlikely(m.size() == 0)) {
29 return ret_type(Eigen::MatrixXd(0, 0));
30 }
31
32 const int M = std::min(m.rows(), m.cols());
33 auto arena_m = to_arena(m);
34
35 Eigen::JacobiSVD<Eigen::MatrixXd> svd(
36 arena_m.val(), Eigen::ComputeThinU | Eigen::ComputeThinV);
37
38 auto arena_D = to_arena(svd.singularValues());
39
40 arena_t<Eigen::MatrixXd> arena_Fm(M, M);
41
42 for (int i = 0; i < M; i++) {
43 for (int j = 0; j < M; j++) {
44 if (j == i) {
45 arena_Fm(i, j) = 0.0;
46 } else {
47 arena_Fm(i, j)
48 = 1.0 / (arena_D[j] - arena_D[i]) - 1.0 / (arena_D[i] + arena_D[j]);
49 }
50 }
51 }
52
53 auto arena_U = to_arena(svd.matrixU());
54 arena_t<ret_type> arena_V = svd.matrixV();
55
56 reverse_pass_callback([arena_m, arena_U, arena_D, arena_V,
57 arena_Fm]() mutable {
58 Eigen::MatrixXd VTVadj = arena_V.val_op().transpose() * arena_V.adj_op();
59 arena_m.adj()
60 += 0.5 * arena_U
61 * (arena_Fm.array() * (VTVadj - VTVadj.transpose()).array())
62 .matrix()
63 * arena_V.val_op().transpose()
64 + arena_U * arena_D.asDiagonal().inverse()
65 * arena_V.adj_op().transpose()
66 * (Eigen::MatrixXd::Identity(arena_m.cols(), arena_m.cols())
67 - arena_V.val_op() * arena_V.val_op().transpose());
68 });
69
70 return ret_type(arena_V);
71}
72
73} // namespace math
74} // namespace stan
75
76#endif
#define unlikely(x)
auto transpose(Arg &&a)
Transposes a kernel generator expression.
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< base_type_t< EigMat >, -1, 1 >, Eigen::Matrix< value_type_t< EigMat >, -1, -1 > > svd(const EigMat &m)
Given input matrix m, return the singular value decomposition (U,D,V) such that m = U*diag(D)*V^{T}
Definition svd.hpp:24
Eigen::Matrix< value_type_t< EigMat >, Eigen::Dynamic, Eigen::Dynamic > svd_V(const EigMat &m)
Given input matrix m, return matrix V where m = UDV^{T}
Definition svd_V.hpp:19
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 ...