Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
SuiteSparseQRSupport.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2012 Desire Nuentsa <[email protected]>
5// Copyright (C) 2014 Gael Guennebaud <[email protected]>
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
17namespace Eigen {
18
19template <typename MatrixType>
20class SPQR;
21template <typename SPQRType>
22struct SPQRMatrixQReturnType;
23template <typename SPQRType>
24struct SPQRMatrixQTransposeReturnType;
25template <typename SPQRType, typename Derived>
26struct SPQR_QProduct;
27namespace internal {
28template <typename SPQRType>
29struct traits<SPQRMatrixQReturnType<SPQRType> > {
30 typedef typename SPQRType::MatrixType ReturnType;
31};
32template <typename SPQRType>
33struct traits<SPQRMatrixQTransposeReturnType<SPQRType> > {
34 typedef typename SPQRType::MatrixType ReturnType;
35};
36template <typename SPQRType, typename Derived>
37struct traits<SPQR_QProduct<SPQRType, Derived> > {
38 typedef typename Derived::PlainObject ReturnType;
39};
40} // End namespace internal
41
66template <typename MatrixType_>
67class SPQR : public SparseSolverBase<SPQR<MatrixType_> > {
68 protected:
70 using Base::m_isInitialized;
71
72 public:
73 typedef typename MatrixType_::Scalar Scalar;
74 typedef typename MatrixType_::RealScalar RealScalar;
75 typedef SuiteSparse_long StorageIndex;
78 enum { ColsAtCompileTime = Dynamic, MaxColsAtCompileTime = Dynamic };
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),
87 m_tolerance(NumTraits<Scalar>::epsilon()),
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),
103 m_tolerance(NumTraits<Scalar>::epsilon()),
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
129 MatrixType mat(matrix);
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) {
149 m_info = NumericalIssue;
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) {
200 m_isRUpToDate = true;
201 }
202 return m_R;
203 }
205 SPQRMatrixQReturnType<SPQR> matrixQ() const { return SPQRMatrixQReturnType<SPQR>(*this); }
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:
241 bool m_analysisIsOk;
242 bool m_factorizationIsOk;
243 mutable bool m_isRUpToDate;
244 mutable ComputationInfo m_info;
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
257 Index m_rows;
258 template <typename, typename>
259 friend struct SPQR_QProduct;
260};
261
262template <typename SPQRType, typename Derived>
263struct 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);
281 res = Matrix<Scalar, ResType::RowsAtCompileTime, ResType::ColsAtCompileTime>::Map(
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;
287 bool m_transpose;
288};
289template <typename SPQRType>
290struct SPQRMatrixQReturnType {
291 SPQRMatrixQReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
292 template <typename Derived>
293 SPQR_QProduct<SPQRType, Derived> operator*(const MatrixBase<Derived>& other) {
294 return SPQR_QProduct<SPQRType, Derived>(m_spqr, other.derived(), false);
295 }
296 SPQRMatrixQTransposeReturnType<SPQRType> adjoint() const { return SPQRMatrixQTransposeReturnType<SPQRType>(m_spqr); }
297 // To use for operations with the transpose of Q
298 SPQRMatrixQTransposeReturnType<SPQRType> transpose() const {
299 return SPQRMatrixQTransposeReturnType<SPQRType>(m_spqr);
300 }
301 const SPQRType& m_spqr;
302};
303
304template <typename SPQRType>
305struct SPQRMatrixQTransposeReturnType {
306 SPQRMatrixQTransposeReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
307 template <typename Derived>
308 SPQR_QProduct<SPQRType, Derived> operator*(const MatrixBase<Derived>& other) {
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
Derived & setZero()
Definition CwiseNullaryOp.h:549
EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition EigenBase.h:61
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
Sparse QR factorization based on SuiteSparseQR library.
Definition SuiteSparseQRSupport.h:67
void setPivotThreshold(const RealScalar &tol)
Set the tolerance tol to treat columns with 2-norm < =tol as zero.
Definition SuiteSparseQRSupport.h:222
Index rank() const
Definition SuiteSparseQRSupport.h:215
cholmod_common * cholmodCommon() const
Definition SuiteSparseQRSupport.h:228
SPQRMatrixQReturnType< SPQR > matrixQ() const
Get an expression of the matrix Q.
Definition SuiteSparseQRSupport.h:205
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SuiteSparseQRSupport.h:235
Index rows() const
Definition SuiteSparseQRSupport.h:160
const MatrixType matrixR() const
Definition SuiteSparseQRSupport.h:196
Index cols() const
Definition SuiteSparseQRSupport.h:165
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
A versatible sparse matrix representation.
Definition SparseUtil.h:47
Index cols() const
Definition SparseMatrix.h:161
Index rows() const
Definition SparseMatrix.h:159
A base class for sparse solvers.
Definition SparseSolverBase.h:67
ComputationInfo
Definition Constants.h:438
@ NumericalIssue
Definition Constants.h:442
@ Success
Definition Constants.h:440
Namespace containing all symbols from the Eigen library.
Definition Core:137
Map< const SparseMatrix< Scalar, ColMajor, StorageIndex > > viewAsEigen(cholmod_sparse &cm)
Definition CholmodSupport.h:153
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
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523