SuiteSparseQRSupport.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Desire Nuentsa <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 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 #ifndef EIGEN_SUITESPARSEQRSUPPORT_H
12 #define EIGEN_SUITESPARSEQRSUPPORT_H
13 
14 // IWYU pragma: private
15 #include "./InternalHeaderCheck.h"
16 
17 namespace Eigen {
18 
19 template <typename MatrixType>
20 class SPQR;
21 template <typename SPQRType>
22 struct SPQRMatrixQReturnType;
23 template <typename SPQRType>
24 struct SPQRMatrixQTransposeReturnType;
25 template <typename SPQRType, typename Derived>
26 struct SPQR_QProduct;
27 namespace internal {
28 template <typename SPQRType>
29 struct traits<SPQRMatrixQReturnType<SPQRType> > {
30  typedef typename SPQRType::MatrixType ReturnType;
31 };
32 template <typename SPQRType>
34  typedef typename SPQRType::MatrixType ReturnType;
35 };
36 template <typename SPQRType, typename Derived>
37 struct traits<SPQR_QProduct<SPQRType, Derived> > {
38  typedef typename Derived::PlainObject ReturnType;
39 };
40 } // End namespace internal
41 
66 template <typename MatrixType_>
67 class SPQR : public SparseSolverBase<SPQR<MatrixType_> > {
68  protected:
71 
72  public:
73  typedef typename MatrixType_::Scalar Scalar;
75  typedef SuiteSparse_long StorageIndex;
79 
80  public:
81  SPQR()
82  : m_analysisIsOk(false),
83  m_factorizationIsOk(false),
84  m_isRUpToDate(false),
85  m_ordering(SPQR_ORDERING_DEFAULT),
86  m_allow_tol(SPQR_DEFAULT_TOL),
88  m_cR(0),
89  m_E(0),
90  m_H(0),
91  m_HPinv(0),
92  m_HTau(0),
93  m_useDefaultThreshold(true) {
94  cholmod_l_start(&m_cc);
95  }
96 
97  explicit SPQR(const MatrixType_& matrix)
98  : m_analysisIsOk(false),
99  m_factorizationIsOk(false),
100  m_isRUpToDate(false),
101  m_ordering(SPQR_ORDERING_DEFAULT),
102  m_allow_tol(SPQR_DEFAULT_TOL),
104  m_cR(0),
105  m_E(0),
106  m_H(0),
107  m_HPinv(0),
108  m_HTau(0),
109  m_useDefaultThreshold(true) {
110  cholmod_l_start(&m_cc);
111  compute(matrix);
112  }
113 
114  ~SPQR() {
115  SPQR_free();
116  cholmod_l_finish(&m_cc);
117  }
118  void SPQR_free() {
119  cholmod_l_free_sparse(&m_H, &m_cc);
120  cholmod_l_free_sparse(&m_cR, &m_cc);
121  cholmod_l_free_dense(&m_HTau, &m_cc);
122  std::free(m_E);
123  std::free(m_HPinv);
124  }
125 
126  void compute(const MatrixType_& matrix) {
127  if (m_isInitialized) SPQR_free();
128 
130 
131  /* Compute the default threshold as in MatLab, see:
132  * Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing
133  * Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011, Page 8:3
134  */
135  RealScalar pivotThreshold = m_tolerance;
136  if (m_useDefaultThreshold) {
137  RealScalar max2Norm = 0.0;
138  for (int j = 0; j < mat.cols(); j++) max2Norm = numext::maxi(max2Norm, mat.col(j).norm());
139  if (numext::is_exactly_zero(max2Norm)) max2Norm = RealScalar(1);
140  pivotThreshold = 20 * (mat.rows() + mat.cols()) * max2Norm * NumTraits<RealScalar>::epsilon();
141  }
142  cholmod_sparse A;
143  A = viewAsCholmod(mat);
144  m_rows = matrix.rows();
145  m_rank = SuiteSparseQR<Scalar>(m_ordering, pivotThreshold, internal::convert_index<StorageIndex>(matrix.cols()), &A,
146  &m_cR, &m_E, &m_H, &m_HPinv, &m_HTau, &m_cc);
147 
148  if (!m_cR) {
150  m_isInitialized = false;
151  return;
152  }
153  m_info = Success;
154  m_isInitialized = true;
155  m_isRUpToDate = false;
156  }
160  inline Index rows() const { return m_rows; }
161 
165  inline Index cols() const { return m_cR->ncol; }
166 
167  template <typename Rhs, typename Dest>
168  void _solve_impl(const MatrixBase<Rhs>& b, MatrixBase<Dest>& dest) const {
169  eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
170  eigen_assert(b.cols() == 1 && "This method is for vectors only");
171 
172  // Compute Q^T * b
173  typename Dest::PlainObject y, y2;
174  y = matrixQ().transpose() * b;
175 
176  // Solves with the triangular matrix R
177  Index rk = this->rank();
178  y2 = y;
179  y.resize((std::max)(cols(), Index(y.rows())), y.cols());
180  y.topRows(rk) = this->matrixR().topLeftCorner(rk, rk).template triangularView<Upper>().solve(y2.topRows(rk));
181 
182  // Apply the column permutation
183  // colsPermutation() performs a copy of the permutation,
184  // so let's apply it manually:
185  for (Index i = 0; i < rk; ++i) dest.row(m_E[i]) = y.row(i);
186  for (Index i = rk; i < cols(); ++i) dest.row(m_E[i]).setZero();
187 
188  // y.bottomRows(y.rows()-rk).setZero();
189  // dest = colsPermutation() * y.topRows(cols());
190 
191  m_info = Success;
192  }
193 
196  const MatrixType matrixR() const {
197  eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
198  if (!m_isRUpToDate) {
199  m_R = viewAsEigen<Scalar, StorageIndex>(*m_cR);
200  m_isRUpToDate = true;
201  }
202  return m_R;
203  }
208  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
209  return PermutationType(m_E, m_cR->ncol);
210  }
215  Index rank() const {
216  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
217  return m_cc.SPQR_istat[4];
218  }
220  void setSPQROrdering(int ord) { m_ordering = ord; }
222  void setPivotThreshold(const RealScalar& tol) {
223  m_useDefaultThreshold = false;
224  m_tolerance = tol;
225  }
226 
228  cholmod_common* cholmodCommon() const { return &m_cc; }
229 
236  eigen_assert(m_isInitialized && "Decomposition is not initialized.");
237  return m_info;
238  }
239 
240  protected:
243  mutable bool m_isRUpToDate;
245  int m_ordering; // Ordering method to use, see SPQR's manual
246  int m_allow_tol; // Allow to use some tolerance during numerical factorization.
247  RealScalar m_tolerance; // treat columns with 2-norm below this tolerance as zero
248  mutable cholmod_sparse* m_cR = nullptr; // The sparse R factor in cholmod format
249  mutable MatrixType m_R; // The sparse matrix R in Eigen format
250  mutable StorageIndex* m_E = nullptr; // The permutation applied to columns
251  mutable cholmod_sparse* m_H = nullptr; // The householder vectors
252  mutable StorageIndex* m_HPinv = nullptr; // The row permutation of H
253  mutable cholmod_dense* m_HTau = nullptr; // The Householder coefficients
254  mutable Index m_rank; // The rank of the matrix
255  mutable cholmod_common m_cc; // Workspace and parameters
256  bool m_useDefaultThreshold; // Use default threshold
258  template <typename, typename>
259  friend struct SPQR_QProduct;
260 };
261 
262 template <typename SPQRType, typename Derived>
263 struct SPQR_QProduct : ReturnByValue<SPQR_QProduct<SPQRType, Derived> > {
264  typedef typename SPQRType::Scalar Scalar;
265  typedef typename SPQRType::StorageIndex StorageIndex;
266  // Define the constructor to get reference to argument types
267  SPQR_QProduct(const SPQRType& spqr, const Derived& other, bool transpose)
268  : m_spqr(spqr), m_other(other), m_transpose(transpose) {}
269 
270  inline Index rows() const { return m_transpose ? m_spqr.rows() : m_spqr.cols(); }
271  inline Index cols() const { return m_other.cols(); }
272  // Assign to a vector
273  template <typename ResType>
274  void evalTo(ResType& res) const {
275  cholmod_dense y_cd;
276  cholmod_dense* x_cd;
277  int method = m_transpose ? SPQR_QTX : SPQR_QX;
278  cholmod_common* cc = m_spqr.cholmodCommon();
279  y_cd = viewAsCholmod(m_other.const_cast_derived());
280  x_cd = SuiteSparseQR_qmult<Scalar>(method, m_spqr.m_H, m_spqr.m_HTau, m_spqr.m_HPinv, &y_cd, cc);
282  reinterpret_cast<Scalar*>(x_cd->x), x_cd->nrow, x_cd->ncol);
283  cholmod_l_free_dense(&x_cd, cc);
284  }
285  const SPQRType& m_spqr;
286  const Derived& m_other;
288 };
289 template <typename SPQRType>
291  SPQRMatrixQReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
292  template <typename Derived>
294  return SPQR_QProduct<SPQRType, Derived>(m_spqr, other.derived(), false);
295  }
297  // To use for operations with the transpose of Q
300  }
301  const SPQRType& m_spqr;
302 };
303 
304 template <typename SPQRType>
306  SPQRMatrixQTransposeReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
307  template <typename Derived>
309  return SPQR_QProduct<SPQRType, Derived>(m_spqr, other.derived(), true);
310  }
311  const SPQRType& m_spqr;
312 };
313 
314 } // End namespace Eigen
315 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
Eigen::SparseMatrix< double > mat
Definition: EigenUnitTest.cpp:10
#define eigen_assert(x)
Definition: Macros.h:910
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition: PartialRedux_count.cpp:3
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
EIGEN_DEVICE_FUNC Derived & setZero()
Definition: CwiseNullaryOp.h:554
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
static ConstMapType Map(const Scalar *data)
Definition: PlainObjectBase.h:595
Definition: ReturnByValue.h:50
Sparse QR factorization based on SuiteSparseQR library.
Definition: SuiteSparseQRSupport.h:67
cholmod_dense * m_HTau
Definition: SuiteSparseQRSupport.h:253
bool m_factorizationIsOk
Definition: SuiteSparseQRSupport.h:242
void compute(const MatrixType_ &matrix)
Definition: SuiteSparseQRSupport.h:126
void setPivotThreshold(const RealScalar &tol)
Set the tolerance tol to treat columns with 2-norm < =tol as zero.
Definition: SuiteSparseQRSupport.h:222
StorageIndex * m_E
Definition: SuiteSparseQRSupport.h:250
void _solve_impl(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
Definition: SuiteSparseQRSupport.h:168
SparseMatrix< Scalar, ColMajor, StorageIndex > MatrixType
Definition: SuiteSparseQRSupport.h:76
SuiteSparse_long StorageIndex
Definition: SuiteSparseQRSupport.h:75
cholmod_sparse * m_cR
Definition: SuiteSparseQRSupport.h:248
SPQR(const MatrixType_ &matrix)
Definition: SuiteSparseQRSupport.h:97
Index rank() const
Definition: SuiteSparseQRSupport.h:215
~SPQR()
Definition: SuiteSparseQRSupport.h:114
bool m_isRUpToDate
Definition: SuiteSparseQRSupport.h:243
SPQR()
Definition: SuiteSparseQRSupport.h:81
cholmod_common * cholmodCommon() const
Definition: SuiteSparseQRSupport.h:228
ComputationInfo info() const
Reports whether previous computation was successful.
Definition: SuiteSparseQRSupport.h:235
Index m_rows
Definition: SuiteSparseQRSupport.h:257
RealScalar m_tolerance
Definition: SuiteSparseQRSupport.h:247
Index rows() const
Definition: SuiteSparseQRSupport.h:160
MatrixType_::RealScalar RealScalar
Definition: SuiteSparseQRSupport.h:74
int m_ordering
Definition: SuiteSparseQRSupport.h:245
cholmod_sparse * m_H
Definition: SuiteSparseQRSupport.h:251
const MatrixType matrixR() const
Definition: SuiteSparseQRSupport.h:196
MatrixType_::Scalar Scalar
Definition: SuiteSparseQRSupport.h:73
Index m_rank
Definition: SuiteSparseQRSupport.h:254
bool m_useDefaultThreshold
Definition: SuiteSparseQRSupport.h:256
void SPQR_free()
Definition: SuiteSparseQRSupport.h:118
StorageIndex * m_HPinv
Definition: SuiteSparseQRSupport.h:252
int m_allow_tol
Definition: SuiteSparseQRSupport.h:246
ComputationInfo m_info
Definition: SuiteSparseQRSupport.h:244
Index cols() const
Definition: SuiteSparseQRSupport.h:165
Map< PermutationMatrix< Dynamic, Dynamic, StorageIndex > > PermutationType
Definition: SuiteSparseQRSupport.h:77
bool m_analysisIsOk
Definition: SuiteSparseQRSupport.h:241
MatrixType m_R
Definition: SuiteSparseQRSupport.h:249
SparseSolverBase< SPQR< MatrixType_ > > Base
Definition: SuiteSparseQRSupport.h:69
@ MaxColsAtCompileTime
Definition: SuiteSparseQRSupport.h:78
@ ColsAtCompileTime
Definition: SuiteSparseQRSupport.h:78
cholmod_common m_cc
Definition: SuiteSparseQRSupport.h:255
SPQRMatrixQReturnType< SPQR > matrixQ() const
Get an expression of the matrix Q.
Definition: SuiteSparseQRSupport.h:205
PermutationType colsPermutation() const
Get the permutation that was applied to columns of A.
Definition: SuiteSparseQRSupport.h:207
void setSPQROrdering(int ord)
Set the fill-reducing ordering method to be used.
Definition: SuiteSparseQRSupport.h:220
RealScalar norm() const
Definition: SparseDot.h:88
Index cols() const
Definition: SparseMatrix.h:161
Index rows() const
Definition: SparseMatrix.h:159
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
bool m_isInitialized
Definition: SparseSolverBase.h:110
Eigen::Map< Eigen::Matrix< T, Eigen::Dynamic, Eigen::Dynamic, Eigen::ColMajor >, 0, Eigen::OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition: common.h:85
#define max(a, b)
Definition: datatypes.h:23
ComputationInfo
Definition: Constants.h:438
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
Scalar * y
Definition: level1_cplx_impl.h:128
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition: MathFunctions.h:926
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC bool is_exactly_zero(const X &x)
Definition: Meta.h:592
Namespace containing all symbols from the Eigen library.
Definition: bench_norm.cpp:70
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< Scalar_, Options_, StorageIndex_ > > mat)
Definition: CholmodSupport.h:64
const int Dynamic
Definition: Constants.h:25
void transpose()
Definition: skew_symmetric_matrix3.cpp:135
Definition: Eigen_Colamd.h:49
double epsilon
Definition: osc_ring_sarah_asymptotics.h:43
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: SuiteSparseQRSupport.h:290
SPQRMatrixQTransposeReturnType< SPQRType > adjoint() const
Definition: SuiteSparseQRSupport.h:296
SPQRMatrixQReturnType(const SPQRType &spqr)
Definition: SuiteSparseQRSupport.h:291
const SPQRType & m_spqr
Definition: SuiteSparseQRSupport.h:301
SPQRMatrixQTransposeReturnType< SPQRType > transpose() const
Definition: SuiteSparseQRSupport.h:298
SPQR_QProduct< SPQRType, Derived > operator*(const MatrixBase< Derived > &other)
Definition: SuiteSparseQRSupport.h:293
Definition: SuiteSparseQRSupport.h:305
SPQRMatrixQTransposeReturnType(const SPQRType &spqr)
Definition: SuiteSparseQRSupport.h:306
const SPQRType & m_spqr
Definition: SuiteSparseQRSupport.h:311
SPQR_QProduct< SPQRType, Derived > operator*(const MatrixBase< Derived > &other)
Definition: SuiteSparseQRSupport.h:308
Definition: SuiteSparseQRSupport.h:263
bool m_transpose
Definition: SuiteSparseQRSupport.h:287
SPQRType::Scalar Scalar
Definition: SuiteSparseQRSupport.h:264
SPQR_QProduct(const SPQRType &spqr, const Derived &other, bool transpose)
Definition: SuiteSparseQRSupport.h:267
void evalTo(ResType &res) const
Definition: SuiteSparseQRSupport.h:274
const SPQRType & m_spqr
Definition: SuiteSparseQRSupport.h:285
SPQRType::StorageIndex StorageIndex
Definition: SuiteSparseQRSupport.h:265
Index cols() const
Definition: SuiteSparseQRSupport.h:271
Index rows() const
Definition: SuiteSparseQRSupport.h:270
const Derived & m_other
Definition: SuiteSparseQRSupport.h:286
SPQRType::MatrixType ReturnType
Definition: SuiteSparseQRSupport.h:30
SPQRType::MatrixType ReturnType
Definition: SuiteSparseQRSupport.h:34
Derived::PlainObject ReturnType
Definition: SuiteSparseQRSupport.h:38
Definition: ForwardDeclarations.h:21
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2