xref: /aosp_15_r20/external/eigen/test/sparseqr.cpp (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2012 Desire Nuentsa Wakam <[email protected]>
5*bf2c3715SXin Li // Copyright (C) 2014 Gael Guennebaud <[email protected]>
6*bf2c3715SXin Li //
7*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
8*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
9*bf2c3715SXin Li #include "sparse.h"
10*bf2c3715SXin Li #include <Eigen/SparseQR>
11*bf2c3715SXin Li 
12*bf2c3715SXin Li template<typename MatrixType,typename DenseMat>
generate_sparse_rectangular_problem(MatrixType & A,DenseMat & dA,int maxRows=300,int maxCols=150)13*bf2c3715SXin Li int generate_sparse_rectangular_problem(MatrixType& A, DenseMat& dA, int maxRows = 300, int maxCols = 150)
14*bf2c3715SXin Li {
15*bf2c3715SXin Li   eigen_assert(maxRows >= maxCols);
16*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
17*bf2c3715SXin Li   int rows = internal::random<int>(1,maxRows);
18*bf2c3715SXin Li   int cols = internal::random<int>(1,maxCols);
19*bf2c3715SXin Li   double density = (std::max)(8./(rows*cols), 0.01);
20*bf2c3715SXin Li 
21*bf2c3715SXin Li   A.resize(rows,cols);
22*bf2c3715SXin Li   dA.resize(rows,cols);
23*bf2c3715SXin Li   initSparse<Scalar>(density, dA, A,ForceNonZeroDiag);
24*bf2c3715SXin Li   A.makeCompressed();
25*bf2c3715SXin Li   int nop = internal::random<int>(0, internal::random<double>(0,1) > 0.5 ? cols/2 : 0);
26*bf2c3715SXin Li   for(int k=0; k<nop; ++k)
27*bf2c3715SXin Li   {
28*bf2c3715SXin Li     int j0 = internal::random<int>(0,cols-1);
29*bf2c3715SXin Li     int j1 = internal::random<int>(0,cols-1);
30*bf2c3715SXin Li     Scalar s = internal::random<Scalar>();
31*bf2c3715SXin Li     A.col(j0)  = s * A.col(j1);
32*bf2c3715SXin Li     dA.col(j0) = s * dA.col(j1);
33*bf2c3715SXin Li   }
34*bf2c3715SXin Li 
35*bf2c3715SXin Li //   if(rows<cols) {
36*bf2c3715SXin Li //     A.conservativeResize(cols,cols);
37*bf2c3715SXin Li //     dA.conservativeResize(cols,cols);
38*bf2c3715SXin Li //     dA.bottomRows(cols-rows).setZero();
39*bf2c3715SXin Li //   }
40*bf2c3715SXin Li 
41*bf2c3715SXin Li   return rows;
42*bf2c3715SXin Li }
43*bf2c3715SXin Li 
test_sparseqr_scalar()44*bf2c3715SXin Li template<typename Scalar> void test_sparseqr_scalar()
45*bf2c3715SXin Li {
46*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
47*bf2c3715SXin Li   typedef SparseMatrix<Scalar,ColMajor> MatrixType;
48*bf2c3715SXin Li   typedef Matrix<Scalar,Dynamic,Dynamic> DenseMat;
49*bf2c3715SXin Li   typedef Matrix<Scalar,Dynamic,1> DenseVector;
50*bf2c3715SXin Li   MatrixType A;
51*bf2c3715SXin Li   DenseMat dA;
52*bf2c3715SXin Li   DenseVector refX,x,b;
53*bf2c3715SXin Li   SparseQR<MatrixType, COLAMDOrdering<int> > solver;
54*bf2c3715SXin Li   generate_sparse_rectangular_problem(A,dA);
55*bf2c3715SXin Li 
56*bf2c3715SXin Li   b = dA * DenseVector::Random(A.cols());
57*bf2c3715SXin Li   solver.compute(A);
58*bf2c3715SXin Li 
59*bf2c3715SXin Li   // Q should be MxM
60*bf2c3715SXin Li   VERIFY_IS_EQUAL(solver.matrixQ().rows(), A.rows());
61*bf2c3715SXin Li   VERIFY_IS_EQUAL(solver.matrixQ().cols(), A.rows());
62*bf2c3715SXin Li 
63*bf2c3715SXin Li   // R should be MxN
64*bf2c3715SXin Li   VERIFY_IS_EQUAL(solver.matrixR().rows(), A.rows());
65*bf2c3715SXin Li   VERIFY_IS_EQUAL(solver.matrixR().cols(), A.cols());
66*bf2c3715SXin Li 
67*bf2c3715SXin Li   // Q and R can be multiplied
68*bf2c3715SXin Li   DenseMat recoveredA = solver.matrixQ()
69*bf2c3715SXin Li                       * DenseMat(solver.matrixR().template triangularView<Upper>())
70*bf2c3715SXin Li                       * solver.colsPermutation().transpose();
71*bf2c3715SXin Li   VERIFY_IS_EQUAL(recoveredA.rows(), A.rows());
72*bf2c3715SXin Li   VERIFY_IS_EQUAL(recoveredA.cols(), A.cols());
73*bf2c3715SXin Li 
74*bf2c3715SXin Li   // and in the full rank case the original matrix is recovered
75*bf2c3715SXin Li   if (solver.rank() == A.cols())
76*bf2c3715SXin Li   {
77*bf2c3715SXin Li       VERIFY_IS_APPROX(A, recoveredA);
78*bf2c3715SXin Li   }
79*bf2c3715SXin Li 
80*bf2c3715SXin Li   if(internal::random<float>(0,1)>0.5f)
81*bf2c3715SXin Li     solver.factorize(A);  // this checks that calling analyzePattern is not needed if the pattern do not change.
82*bf2c3715SXin Li   if (solver.info() != Success)
83*bf2c3715SXin Li   {
84*bf2c3715SXin Li     std::cerr << "sparse QR factorization failed\n";
85*bf2c3715SXin Li     exit(0);
86*bf2c3715SXin Li     return;
87*bf2c3715SXin Li   }
88*bf2c3715SXin Li   x = solver.solve(b);
89*bf2c3715SXin Li   if (solver.info() != Success)
90*bf2c3715SXin Li   {
91*bf2c3715SXin Li     std::cerr << "sparse QR factorization failed\n";
92*bf2c3715SXin Li     exit(0);
93*bf2c3715SXin Li     return;
94*bf2c3715SXin Li   }
95*bf2c3715SXin Li 
96*bf2c3715SXin Li   // Compare with a dense QR solver
97*bf2c3715SXin Li   ColPivHouseholderQR<DenseMat> dqr(dA);
98*bf2c3715SXin Li   refX = dqr.solve(b);
99*bf2c3715SXin Li 
100*bf2c3715SXin Li   bool rank_deficient = A.cols()>A.rows() || dqr.rank()<A.cols();
101*bf2c3715SXin Li   if(rank_deficient)
102*bf2c3715SXin Li   {
103*bf2c3715SXin Li     // rank deficient problem -> we might have to increase the threshold
104*bf2c3715SXin Li     // to get a correct solution.
105*bf2c3715SXin Li     RealScalar th = RealScalar(20)*dA.colwise().norm().maxCoeff()*(A.rows()+A.cols()) * NumTraits<RealScalar>::epsilon();
106*bf2c3715SXin Li     for(Index k=0; (k<16) && !test_isApprox(A*x,b); ++k)
107*bf2c3715SXin Li     {
108*bf2c3715SXin Li       th *= RealScalar(10);
109*bf2c3715SXin Li       solver.setPivotThreshold(th);
110*bf2c3715SXin Li       solver.compute(A);
111*bf2c3715SXin Li       x = solver.solve(b);
112*bf2c3715SXin Li     }
113*bf2c3715SXin Li   }
114*bf2c3715SXin Li 
115*bf2c3715SXin Li   VERIFY_IS_APPROX(A * x, b);
116*bf2c3715SXin Li 
117*bf2c3715SXin Li   // For rank deficient problem, the estimated rank might
118*bf2c3715SXin Li   // be slightly off, so let's only raise a warning in such cases.
119*bf2c3715SXin Li   if(rank_deficient) ++g_test_level;
120*bf2c3715SXin Li   VERIFY_IS_EQUAL(solver.rank(), dqr.rank());
121*bf2c3715SXin Li   if(rank_deficient) --g_test_level;
122*bf2c3715SXin Li 
123*bf2c3715SXin Li   if(solver.rank()==A.cols()) // full rank
124*bf2c3715SXin Li     VERIFY_IS_APPROX(x, refX);
125*bf2c3715SXin Li //   else
126*bf2c3715SXin Li //     VERIFY((dA * refX - b).norm() * 2 > (A * x - b).norm() );
127*bf2c3715SXin Li 
128*bf2c3715SXin Li   // Compute explicitly the matrix Q
129*bf2c3715SXin Li   MatrixType Q, QtQ, idM;
130*bf2c3715SXin Li   Q = solver.matrixQ();
131*bf2c3715SXin Li   //Check  ||Q' * Q - I ||
132*bf2c3715SXin Li   QtQ = Q * Q.adjoint();
133*bf2c3715SXin Li   idM.resize(Q.rows(), Q.rows()); idM.setIdentity();
134*bf2c3715SXin Li   VERIFY(idM.isApprox(QtQ));
135*bf2c3715SXin Li 
136*bf2c3715SXin Li   // Q to dense
137*bf2c3715SXin Li   DenseMat dQ;
138*bf2c3715SXin Li   dQ = solver.matrixQ();
139*bf2c3715SXin Li   VERIFY_IS_APPROX(Q, dQ);
140*bf2c3715SXin Li }
EIGEN_DECLARE_TEST(sparseqr)141*bf2c3715SXin Li EIGEN_DECLARE_TEST(sparseqr)
142*bf2c3715SXin Li {
143*bf2c3715SXin Li   for(int i=0; i<g_repeat; ++i)
144*bf2c3715SXin Li   {
145*bf2c3715SXin Li     CALL_SUBTEST_1(test_sparseqr_scalar<double>());
146*bf2c3715SXin Li     CALL_SUBTEST_2(test_sparseqr_scalar<std::complex<double> >());
147*bf2c3715SXin Li   }
148*bf2c3715SXin Li }
149*bf2c3715SXin Li 
150