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