Automatic Differentiation
 
Loading...
Searching...
No Matches
generalized_inverse.hpp
Go to the documentation of this file.
1#ifndef STAN_MATH_REV_FUN_GENERALIZED_INVERSE_HPP
2#define STAN_MATH_REV_FUN_GENERALIZED_INVERSE_HPP
3
10
11namespace stan {
12namespace math {
13
14namespace internal {
15/*
16 * Reverse mode specialization of calculating the generalized inverse of a
17 * matrix.
18 * <ul><li> Golub, G.H. and Pereyra, V. The Differentiation of Pseudo-Inverses
19 * and Nonlinear Least Squares Problems Whose Variables Separate. <i>SIAM
20 * Journal on Numerical Analysis</i>, Vol. 10, No. 2 (Apr., 1973), pp.
21 * 413-432</li></ul>
22 */
23template <typename T1, typename T2>
24inline auto generalized_inverse_lambda(T1& G_arena, T2& inv_G) {
25 return [G_arena, inv_G]() mutable {
26 G_arena.adj()
27 += -(inv_G.val_op().transpose() * inv_G.adj_op()
28 * inv_G.val_op().transpose())
29 + (-G_arena.val_op() * inv_G.val_op()
30 + Eigen::MatrixXd::Identity(G_arena.rows(), inv_G.cols()))
31 * inv_G.adj_op().transpose() * inv_G.val_op()
32 * inv_G.val_op().transpose()
33 + inv_G.val_op().transpose() * inv_G.val_op()
34 * inv_G.adj_op().transpose()
35 * (-inv_G.val_op() * G_arena.val_op()
36 + Eigen::MatrixXd::Identity(inv_G.rows(), G_arena.cols()));
37 };
38}
39} // namespace internal
40
41/*
42 * Reverse mode specialization of calculating the generalized inverse of a
43 * matrix.
44 *
45 * @param G specified matrix
46 * @return Generalized inverse of the matrix (an empty matrix if the specified
47 * matrix has size zero).
48 *
49 * @note For the derivatives of this function to exist the matrix must be
50 * of constant rank.
51 * Reverse mode differentiation algorithm reference:
52 *
53 * <ul><li> Golub, G.H. and Pereyra, V. The Differentiation of Pseudo-Inverses
54 * and Nonlinear Least Squares Problems Whose Variables Separate. <i>SIAM
55 * Journal on Numerical Analysis</i>, Vol. 10, No. 2 (Apr., 1973), pp.
56 * 413-432</li></ul>
57 *
58 * Equation 4.12 in the paper
59 *
60 * See also
61 * http://mathoverflow.net/questions/25778/analytical-formula-for-numerical-derivative-of-the-matrix-pseudo-inverse
62 *
63 */
64template <typename VarMat, require_rev_matrix_t<VarMat>* = nullptr>
65inline auto generalized_inverse(const VarMat& G) {
67
68 if (G.size() == 0)
69 return ret_type(G);
70
71 if (G.rows() == G.cols()) {
72 arena_t<VarMat> G_arena(G);
73 Eigen::CompleteOrthogonalDecomposition<Eigen::MatrixXd>
74 complete_ortho_decomp_G
75 = G_arena.val().completeOrthogonalDecomposition();
76 if (!(complete_ortho_decomp_G.rank() < G.rows())) {
77 return ret_type(inverse(G_arena));
78 } else {
79 arena_t<ret_type> inv_G(complete_ortho_decomp_G.pseudoInverse());
82 return ret_type(inv_G);
83 }
84 } else if (G.rows() < G.cols()) {
85 arena_t<VarMat> G_arena(G);
86 arena_t<ret_type> inv_G((G_arena.val_op() * G_arena.val_op().transpose())
87 .ldlt()
88 .solve(G_arena.val_op())
89 .transpose());
91 return ret_type(inv_G);
92 } else {
93 arena_t<VarMat> G_arena(G);
94 arena_t<ret_type> inv_G((G_arena.val_op().transpose() * G_arena.val_op())
95 .ldlt()
96 .solve(G_arena.val_op().transpose()));
98 return ret_type(inv_G);
99 }
100}
101
102} // namespace math
103} // namespace stan
104#endif
auto generalized_inverse_lambda(T1 &G_arena, T2 &inv_G)
void reverse_pass_callback(F &&functor)
Puts a callback on the autodiff stack to be called in reverse pass.
Eigen::Matrix< value_type_t< EigMat >, EigMat::ColsAtCompileTime, EigMat::RowsAtCompileTime > generalized_inverse(const EigMat &G)
Returns the Moore-Penrose generalized inverse of the specified matrix.
Eigen::Matrix< value_type_t< EigMat >, EigMat::RowsAtCompileTime, EigMat::ColsAtCompileTime > inverse(const EigMat &m)
Forward mode specialization of calculating the inverse of the matrix.
Definition inverse.hpp:28
std::conditional_t< is_any_var_matrix< ReturnType, Types... >::value, stan::math::var_value< stan::math::promote_scalar_t< double, ReturnType > >, stan::math::promote_scalar_t< stan::math::var_value< double >, ReturnType > > promote_var_matrix_t
Given an Eigen type and several inputs, determine if a matrix should be var<Matrix> or Matrix<var>.
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.
The lgamma implementation in stan-math is based on either the reentrant safe lgamma_r implementation ...