blob: f0e721fcea2cc6949ab8b4ac66a0f3434b6c7dc5 [file] [log] [blame]
Brian Silverman72890c22015-09-19 14:37:37 -04001// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2012 Desire Nuentsa Wakam <desire.nuentsa_wakam@inria.fr>
5// Copyright (C) 2014 Gael Guennebaud <gael.guennebaud@inria.fr>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9#include "sparse.h"
10#include <Eigen/SparseQR>
11
12template<typename MatrixType,typename DenseMat>
Austin Schuh189376f2018-12-20 22:11:15 +110013int generate_sparse_rectangular_problem(MatrixType& A, DenseMat& dA, int maxRows = 300, int maxCols = 150)
Brian Silverman72890c22015-09-19 14:37:37 -040014{
Austin Schuh189376f2018-12-20 22:11:15 +110015 eigen_assert(maxRows >= maxCols);
Brian Silverman72890c22015-09-19 14:37:37 -040016 typedef typename MatrixType::Scalar Scalar;
17 int rows = internal::random<int>(1,maxRows);
Austin Schuh189376f2018-12-20 22:11:15 +110018 int cols = internal::random<int>(1,maxCols);
Brian Silverman72890c22015-09-19 14:37:37 -040019 double density = (std::max)(8./(rows*cols), 0.01);
20
21 A.resize(rows,cols);
22 dA.resize(rows,cols);
23 initSparse<Scalar>(density, dA, A,ForceNonZeroDiag);
24 A.makeCompressed();
25 int nop = internal::random<int>(0, internal::random<double>(0,1) > 0.5 ? cols/2 : 0);
26 for(int k=0; k<nop; ++k)
27 {
28 int j0 = internal::random<int>(0,cols-1);
29 int j1 = internal::random<int>(0,cols-1);
30 Scalar s = internal::random<Scalar>();
31 A.col(j0) = s * A.col(j1);
32 dA.col(j0) = s * dA.col(j1);
33 }
34
35// if(rows<cols) {
36// A.conservativeResize(cols,cols);
37// dA.conservativeResize(cols,cols);
38// dA.bottomRows(cols-rows).setZero();
39// }
40
41 return rows;
42}
43
44template<typename Scalar> void test_sparseqr_scalar()
45{
46 typedef SparseMatrix<Scalar,ColMajor> MatrixType;
47 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMat;
48 typedef Matrix<Scalar,Dynamic,1> DenseVector;
49 MatrixType A;
50 DenseMat dA;
51 DenseVector refX,x,b;
52 SparseQR<MatrixType, COLAMDOrdering<int> > solver;
53 generate_sparse_rectangular_problem(A,dA);
54
55 b = dA * DenseVector::Random(A.cols());
56 solver.compute(A);
Austin Schuh189376f2018-12-20 22:11:15 +110057
58 // Q should be MxM
59 VERIFY_IS_EQUAL(solver.matrixQ().rows(), A.rows());
60 VERIFY_IS_EQUAL(solver.matrixQ().cols(), A.rows());
61
62 // R should be MxN
63 VERIFY_IS_EQUAL(solver.matrixR().rows(), A.rows());
64 VERIFY_IS_EQUAL(solver.matrixR().cols(), A.cols());
65
66 // Q and R can be multiplied
67 DenseMat recoveredA = solver.matrixQ()
68 * DenseMat(solver.matrixR().template triangularView<Upper>())
69 * solver.colsPermutation().transpose();
70 VERIFY_IS_EQUAL(recoveredA.rows(), A.rows());
71 VERIFY_IS_EQUAL(recoveredA.cols(), A.cols());
72
73 // and in the full rank case the original matrix is recovered
74 if (solver.rank() == A.cols())
75 {
76 VERIFY_IS_APPROX(A, recoveredA);
77 }
78
79 if(internal::random<float>(0,1)>0.5f)
Brian Silverman72890c22015-09-19 14:37:37 -040080 solver.factorize(A); // this checks that calling analyzePattern is not needed if the pattern do not change.
81 if (solver.info() != Success)
82 {
83 std::cerr << "sparse QR factorization failed\n";
84 exit(0);
85 return;
86 }
87 x = solver.solve(b);
88 if (solver.info() != Success)
89 {
90 std::cerr << "sparse QR factorization failed\n";
91 exit(0);
92 return;
93 }
94
95 VERIFY_IS_APPROX(A * x, b);
96
97 //Compare with a dense QR solver
98 ColPivHouseholderQR<DenseMat> dqr(dA);
99 refX = dqr.solve(b);
100
101 VERIFY_IS_EQUAL(dqr.rank(), solver.rank());
102 if(solver.rank()==A.cols()) // full rank
103 VERIFY_IS_APPROX(x, refX);
104// else
105// VERIFY((dA * refX - b).norm() * 2 > (A * x - b).norm() );
106
107 // Compute explicitly the matrix Q
108 MatrixType Q, QtQ, idM;
109 Q = solver.matrixQ();
110 //Check ||Q' * Q - I ||
111 QtQ = Q * Q.adjoint();
112 idM.resize(Q.rows(), Q.rows()); idM.setIdentity();
113 VERIFY(idM.isApprox(QtQ));
Austin Schuh189376f2018-12-20 22:11:15 +1100114
115 // Q to dense
116 DenseMat dQ;
117 dQ = solver.matrixQ();
118 VERIFY_IS_APPROX(Q, dQ);
Brian Silverman72890c22015-09-19 14:37:37 -0400119}
120void test_sparseqr()
121{
122 for(int i=0; i<g_repeat; ++i)
123 {
124 CALL_SUBTEST_1(test_sparseqr_scalar<double>());
125 CALL_SUBTEST_2(test_sparseqr_scalar<std::complex<double> >());
126 }
127}
128