1#ifndef STAN_MATH_REV_FUN_SVD_U_HPP
2#define STAN_MATH_REV_FUN_SVD_U_HPP
25template <
typename EigMat, require_rev_matrix_t<EigMat>* =
nullptr>
26inline auto svd_U(
const EigMat& m) {
29 return ret_type(Eigen::MatrixXd(0, 0));
32 const int M = std::min(m.rows(), m.cols());
35 Eigen::JacobiSVD<Eigen::MatrixXd>
svd(
36 arena_m.val(), Eigen::ComputeThinU | Eigen::ComputeThinV);
42 for (
int i = 0; i < M; i++) {
43 for (
int j = 0; j < M; j++) {
48 = 1.0 / (arena_D[j] - arena_D[i]) + 1.0 / (arena_D[i] + arena_D[j]);
58 Eigen::MatrixXd UUadjT = arena_U.val_op().
transpose() * arena_U.adj_op();
60 += .5 * arena_U.val_op()
61 * (arena_Fp.array() * (UUadjT - UUadjT.transpose()).array())
64 + (Eigen::MatrixXd::Identity(arena_m.rows(), arena_m.rows())
65 - arena_U.val_op() * arena_U.val_op().transpose())
66 * arena_U.adj_op() * arena_D.asDiagonal().inverse()
67 * arena_V.transpose();
70 return ret_type(arena_U);
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}
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...
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}
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 ...