Automatic Differentiation
 
Loading...
Searching...
No Matches
singular_values.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUN_SINGULAR_VALUES_HPP
2#define STAN_MATH_REV_FUN_SINGULAR_VALUES_HPP
3
8
9namespace stan {
10namespace math {
11
22template <typename EigMat, require_rev_matrix_t<EigMat>* = nullptr>
23inline auto singular_values(const EigMat& m) {
25 if (unlikely(m.size() == 0)) {
26 return ret_type(Eigen::VectorXd(0));
27 }
28
29 auto arena_m = to_arena(m);
30
31 Eigen::JacobiSVD<Eigen::MatrixXd> svd(
32 arena_m.val(), Eigen::ComputeThinU | Eigen::ComputeThinV);
33
34 arena_t<ret_type> singular_values = svd.singularValues();
35
36 auto arena_U = to_arena(svd.matrixU());
37 auto arena_V = to_arena(svd.matrixV());
38
39 reverse_pass_callback([arena_m, arena_U, singular_values, arena_V]() mutable {
40 arena_m.adj()
41 += arena_U * singular_values.adj().asDiagonal() * arena_V.transpose();
42 });
43
44 return ret_type(singular_values);
45}
46
47} // namespace math
48} // namespace stan
49
50#endif
#define unlikely(x)
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
auto singular_values(const EigMat &m)
Return the vector of the singular values of the specified matrix in decreasing order of magnitude.
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 ...