Skip to content

Commit d932ecf

Browse files
authored
Merge pull request #2225 from spinkney/feature/ginv
Adding Generalized Inverse
2 parents a87459f + 20742a8 commit d932ecf

File tree

6 files changed

+250
-0
lines changed

6 files changed

+250
-0
lines changed

stan/math/prim/fun.hpp

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -102,6 +102,7 @@
102102
#include <stan/math/prim/fun/fmod.hpp>
103103
#include <stan/math/prim/fun/gamma_p.hpp>
104104
#include <stan/math/prim/fun/gamma_q.hpp>
105+
#include <stan/math/prim/fun/generalized_inverse.hpp>
105106
#include <stan/math/prim/fun/get.hpp>
106107
#include <stan/math/prim/fun/get_base1.hpp>
107108
#include <stan/math/prim/fun/get_base1_lhs.hpp>
Lines changed: 54 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,54 @@
1+
#ifndef STAN_MATH_PRIM_FUN_GENERALIZED_INVERSE_HPP
2+
#define STAN_MATH_PRIM_FUN_GENERALIZED_INVERSE_HPP
3+
4+
#include <stan/math/prim/meta.hpp>
5+
#include <stan/math/prim/err.hpp>
6+
#include <stan/math/prim/fun/Eigen.hpp>
7+
#include <stan/math/prim/fun/to_ref.hpp>
8+
#include <stan/math/prim/fun/inverse.hpp>
9+
10+
namespace stan {
11+
namespace math {
12+
13+
/**
14+
* Returns the Moore-Penrose generalized inverse of the specified matrix.
15+
*
16+
* The method is based on the Cholesky computation of the transform as specified
17+
* in
18+
*
19+
* <ul><li> Courrieu, Pierre. 2008. Fast Computation of Moore-Penrose Inverse
20+
Matrices.
21+
* <i>arXiv</i> <b>0804.4809</b> </li></ul>
22+
*
23+
* @tparam EigMat type of the matrix (must be derived from `Eigen::MatrixBase`)
24+
*
25+
* @param G specified matrix
26+
* @return Generalized inverse of the matrix (an empty matrix if the specified
27+
* matrix has size zero).
28+
*/
29+
template <typename EigMat, require_eigen_t<EigMat>* = nullptr,
30+
require_not_vt_var<EigMat>* = nullptr>
31+
inline Eigen::Matrix<value_type_t<EigMat>, EigMat::ColsAtCompileTime,
32+
EigMat::RowsAtCompileTime>
33+
generalized_inverse(const EigMat& G) {
34+
using value_t = value_type_t<EigMat>;
35+
36+
if (G.size() == 0)
37+
return {};
38+
39+
if (G.rows() == G.cols())
40+
return inverse(G);
41+
42+
const auto& G_ref = to_ref(G);
43+
44+
if (G.rows() < G.cols()) {
45+
return (G_ref * G_ref.transpose()).ldlt().solve(G_ref).transpose();
46+
} else {
47+
return (G_ref.transpose() * G_ref).ldlt().solve(G_ref.transpose());
48+
}
49+
}
50+
51+
} // namespace math
52+
} // namespace stan
53+
54+
#endif

stan/math/rev/fun.hpp

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -62,6 +62,7 @@
6262
#include <stan/math/rev/fun/from_var_value.hpp>
6363
#include <stan/math/rev/fun/gamma_p.hpp>
6464
#include <stan/math/rev/fun/gamma_q.hpp>
65+
#include <stan/math/rev/fun/generalized_inverse.hpp>
6566
#include <stan/math/rev/fun/gp_periodic_cov.hpp>
6667
#include <stan/math/rev/fun/grad.hpp>
6768
#include <stan/math/rev/fun/grad_inc_beta.hpp>
Lines changed: 93 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,93 @@
1+
#ifndef STAN_MATH_REV_FUN_GENERALIZED_INVERSE_HPP
2+
#define STAN_MATH_REV_FUN_GENERALIZED_INVERSE_HPP
3+
4+
#include <stan/math/rev/core.hpp>
5+
#include <stan/math/prim/err.hpp>
6+
#include <stan/math/prim/fun/Eigen.hpp>
7+
#include <stan/math/rev/fun/inverse.hpp>
8+
9+
namespace stan {
10+
namespace math {
11+
12+
namespace internal {
13+
/*
14+
* Reverse mode specialization of calculating the generalized inverse of a
15+
* matrix.
16+
* <ul><li> Golub, G.H. and Pereyra, V. The Differentiation of Pseudo-Inverses
17+
* and Nonlinear Least Squares Problems Whose Variables Separate. <i>SIAM
18+
* Journal on Numerical Analysis</i>, Vol. 10, No. 2 (Apr., 1973), pp.
19+
* 413-432</li></ul>
20+
*/
21+
template <typename T1, typename T2>
22+
inline auto generalized_inverse_lambda(T1& G_arena, T2& inv_G) {
23+
return [G_arena, inv_G]() mutable {
24+
G_arena.adj()
25+
+= -(inv_G.val_op().transpose() * inv_G.adj_op()
26+
* inv_G.val_op().transpose())
27+
+ (-G_arena.val_op() * inv_G.val_op()
28+
+ Eigen::MatrixXd::Identity(G_arena.rows(), inv_G.cols()))
29+
* inv_G.adj_op().transpose() * inv_G.val_op()
30+
* inv_G.val_op().transpose()
31+
+ inv_G.val_op().transpose() * inv_G.val_op()
32+
* inv_G.adj_op().transpose()
33+
* (-inv_G.val_op() * G_arena.val_op()
34+
+ Eigen::MatrixXd::Identity(inv_G.rows(), G_arena.cols()));
35+
};
36+
}
37+
} // namespace internal
38+
39+
/*
40+
* Reverse mode specialization of calculating the generalized inverse of a
41+
* matrix.
42+
*
43+
* @param G specified matrix
44+
* @return Generalized inverse of the matrix (an empty matrix if the specified
45+
* matrix has size zero).
46+
*
47+
* @note For the derivatives of this function to exist the matrix must be
48+
* of constant rank.
49+
* Reverse mode differentiation algorithm reference:
50+
*
51+
* <ul><li> Golub, G.H. and Pereyra, V. The Differentiation of Pseudo-Inverses
52+
* and Nonlinear Least Squares Problems Whose Variables Separate. <i>SIAM
53+
* Journal on Numerical Analysis</i>, Vol. 10, No. 2 (Apr., 1973), pp.
54+
* 413-432</li></ul>
55+
*
56+
* Equation 4.12 in the paper
57+
*
58+
* See also
59+
* http://mathoverflow.net/questions/25778/analytical-formula-for-numerical-derivative-of-the-matrix-pseudo-inverse
60+
*
61+
*/
62+
template <typename VarMat, require_rev_matrix_t<VarMat>* = nullptr>
63+
inline auto generalized_inverse(const VarMat& G) {
64+
using value_t = value_type_t<VarMat>;
65+
using ret_type = promote_var_matrix_t<VarMat, VarMat>;
66+
67+
if (G.size() == 0)
68+
return ret_type(G);
69+
70+
if (G.rows() == G.cols())
71+
return ret_type(inverse(G));
72+
73+
if (G.rows() < G.cols()) {
74+
arena_t<VarMat> G_arena(G);
75+
arena_t<ret_type> inv_G((G_arena.val_op() * G_arena.val_op().transpose())
76+
.ldlt()
77+
.solve(G_arena.val_op())
78+
.transpose());
79+
reverse_pass_callback(internal::generalized_inverse_lambda(G_arena, inv_G));
80+
return ret_type(inv_G);
81+
} else {
82+
arena_t<VarMat> G_arena(G);
83+
arena_t<ret_type> inv_G((G_arena.val_op().transpose() * G_arena.val_op())
84+
.ldlt()
85+
.solve(G_arena.val_op().transpose()));
86+
reverse_pass_callback(internal::generalized_inverse_lambda(G_arena, inv_G));
87+
return ret_type(inv_G);
88+
}
89+
}
90+
91+
} // namespace math
92+
} // namespace stan
93+
#endif
Lines changed: 54 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,54 @@
1+
#include <test/unit/math/test_ad.hpp>
2+
#include <vector>
3+
#include <gtest/gtest.h>
4+
5+
TEST(mathMixMatFun, ad_tests) {
6+
using stan::test::expect_ad;
7+
using stan::test::expect_ad_matvar;
8+
9+
auto f = [](const auto& G) { return stan::math::generalized_inverse(G); };
10+
11+
Eigen::MatrixXd t(0, 0);
12+
expect_ad(f, t);
13+
expect_ad_matvar(f, t);
14+
15+
Eigen::MatrixXd u(1, 1);
16+
u << 2;
17+
expect_ad(f, u);
18+
expect_ad_matvar(f, u);
19+
20+
Eigen::MatrixXd v(2, 3);
21+
v << 1, 3, 5, 2, 4, 6;
22+
expect_ad(f, v);
23+
expect_ad_matvar(f, v);
24+
25+
v << 1.9, 1.3, 2.5, 0.4, 1.7, 0.1;
26+
expect_ad(f, v);
27+
expect_ad_matvar(f, v);
28+
29+
Eigen::MatrixXd s(2, 4);
30+
s << 3.4, 2, 5, 1.2, 2, 1, 3.2, 3.1;
31+
expect_ad(f, s);
32+
expect_ad_matvar(f, s);
33+
34+
// issues around zero require looser tolerances for hessians
35+
stan::test::ad_tolerances tols;
36+
tols.hessian_hessian_ = 2.0;
37+
tols.hessian_fvar_hessian_ = 2.0;
38+
39+
Eigen::MatrixXd w(3, 4);
40+
w << 2, 3, 5, 7, 11, 13, 17, 19, 23, 25, 27, 29;
41+
expect_ad(tols, f, w);
42+
expect_ad_matvar(f, w);
43+
44+
Eigen::MatrixXd z(2, 2);
45+
z << 1, 2, 5, std::numeric_limits<double>::quiet_NaN();
46+
EXPECT_NO_THROW(stan::math::generalized_inverse(z));
47+
48+
// autodiff throws, so following fails (throw behavior must match to pass)
49+
50+
Eigen::MatrixXd a(2, 2);
51+
a << 1.9, 0.3, 0.3, std::numeric_limits<double>::infinity();
52+
expect_ad(f, a);
53+
expect_ad_matvar(f, a);
54+
}
Lines changed: 47 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,47 @@
1+
#include <stan/math/prim.hpp>
2+
#include <test/unit/util.hpp>
3+
#include <gtest/gtest.h>
4+
5+
TEST(MathMatrixPrim, Zero) {
6+
stan::math::matrix_d m0(0, 0);
7+
stan::math::matrix_d ginv = stan::math::generalized_inverse(m0);
8+
EXPECT_EQ(0, ginv.rows());
9+
EXPECT_EQ(0, ginv.cols());
10+
}
11+
12+
TEST(MathMatrixPrim, Singular) {
13+
using stan::math::generalized_inverse;
14+
15+
stan::math::matrix_d m1(2, 3);
16+
m1 << 1, 2, 1, 2, 4, 2;
17+
18+
stan::math::matrix_d m2 = m1 * generalized_inverse(m1) * m1;
19+
EXPECT_MATRIX_NEAR(m1, m2, 1e-9);
20+
}
21+
22+
TEST(MathMatrixPrim, Equal1) {
23+
using stan::math::generalized_inverse;
24+
25+
stan::math::matrix_d m1(2, 3);
26+
m1 << 1, 3, 5, 2, 4, 6;
27+
28+
stan::math::matrix_d m2(2, 2);
29+
m2 << 1, 0, 0, 1;
30+
31+
stan::math::matrix_d m3 = m1 * generalized_inverse(m1);
32+
EXPECT_MATRIX_NEAR(m2, m3, 1e-9);
33+
}
34+
35+
TEST(MathMatrixPrim, Equal2) {
36+
using stan::math::generalized_inverse;
37+
38+
stan::math::matrix_d m1(3, 2);
39+
m1 << 1, 2, 2, 4, 1, 2;
40+
41+
stan::math::matrix_d m2(3, 3);
42+
m2 << 1.0 / 6.0, 1.0 / 3.0, 1.0 / 6.0, 1.0 / 3.0, 2.0 / 3.0, 1.0 / 3.0,
43+
1.0 / 6.0, 1.0 / 3.0, 1.0 / 6.0;
44+
45+
stan::math::matrix_d m3 = m1 * generalized_inverse(m1);
46+
EXPECT_MATRIX_NEAR(m2, m3, 1e-9);
47+
}

0 commit comments

Comments
 (0)