TR-mbed 1.0
Loading...
Searching...
No Matches
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
14namespace Eigen {
15
16 template<typename MatrixType> class SPQR;
17 template<typename SPQRType> struct SPQRMatrixQReturnType;
18 template<typename SPQRType> struct SPQRMatrixQTransposeReturnType;
19 template <typename SPQRType, typename Derived> struct SPQR_QProduct;
20 namespace internal {
21 template <typename SPQRType> struct traits<SPQRMatrixQReturnType<SPQRType> >
22 {
23 typedef typename SPQRType::MatrixType ReturnType;
24 };
25 template <typename SPQRType> struct traits<SPQRMatrixQTransposeReturnType<SPQRType> >
26 {
27 typedef typename SPQRType::MatrixType ReturnType;
28 };
29 template <typename SPQRType, typename Derived> struct traits<SPQR_QProduct<SPQRType, Derived> >
30 {
31 typedef typename Derived::PlainObject ReturnType;
32 };
33 } // End namespace internal
34
59template<typename _MatrixType>
60class SPQR : public SparseSolverBase<SPQR<_MatrixType> >
61{
62 protected:
65 public:
66 typedef typename _MatrixType::Scalar Scalar;
67 typedef typename _MatrixType::RealScalar RealScalar;
68 typedef SuiteSparse_long StorageIndex ;
71 enum {
74 };
75 public:
77 : m_analysisIsOk(false),
79 m_isRUpToDate(false),
80 m_ordering(SPQR_ORDERING_DEFAULT),
81 m_allow_tol(SPQR_DEFAULT_TOL),
82 m_tolerance (NumTraits<Scalar>::epsilon()),
83 m_cR(0),
84 m_E(0),
85 m_H(0),
86 m_HPinv(0),
87 m_HTau(0),
89 {
90 cholmod_l_start(&m_cc);
91 }
92
93 explicit SPQR(const _MatrixType& matrix)
94 : m_analysisIsOk(false),
96 m_isRUpToDate(false),
97 m_ordering(SPQR_ORDERING_DEFAULT),
98 m_allow_tol(SPQR_DEFAULT_TOL),
99 m_tolerance (NumTraits<Scalar>::epsilon()),
100 m_cR(0),
101 m_E(0),
102 m_H(0),
103 m_HPinv(0),
104 m_HTau(0),
106 {
107 cholmod_l_start(&m_cc);
109 }
110
112 {
113 SPQR_free();
114 cholmod_l_finish(&m_cc);
115 }
117 {
118 cholmod_l_free_sparse(&m_H, &m_cc);
119 cholmod_l_free_sparse(&m_cR, &m_cc);
120 cholmod_l_free_dense(&m_HTau, &m_cc);
121 std::free(m_E);
122 std::free(m_HPinv);
123 }
124
125 void compute(const _MatrixType& matrix)
126 {
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;
137 {
138 RealScalar max2Norm = 0.0;
139 for (int j = 0; j < mat.cols(); j++) max2Norm = numext::maxi(max2Norm, mat.col(j).norm());
140 if(max2Norm==RealScalar(0))
141 max2Norm = RealScalar(1);
142 pivotThreshold = 20 * (mat.rows() + mat.cols()) * max2Norm * NumTraits<RealScalar>::epsilon();
143 }
144 cholmod_sparse A;
146 m_rows = matrix.rows();
147 Index col = matrix.cols();
148 m_rank = SuiteSparseQR<Scalar>(m_ordering, pivotThreshold, col, &A,
149 &m_cR, &m_E, &m_H, &m_HPinv, &m_HTau, &m_cc);
150
151 if (!m_cR)
152 {
154 m_isInitialized = false;
155 return;
156 }
157 m_info = Success;
158 m_isInitialized = true;
159 m_isRUpToDate = false;
160 }
164 inline Index rows() const {return m_rows; }
165
169 inline Index cols() const { return m_cR->ncol; }
170
171 template<typename Rhs, typename Dest>
173 {
174 eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
175 eigen_assert(b.cols()==1 && "This method is for vectors only");
176
177 //Compute Q^T * b
178 typename Dest::PlainObject y, y2;
179 y = matrixQ().transpose() * b;
180
181 // Solves with the triangular matrix R
182 Index rk = this->rank();
183 y2 = y;
184 y.resize((std::max)(cols(),Index(y.rows())),y.cols());
185 y.topRows(rk) = this->matrixR().topLeftCorner(rk, rk).template triangularView<Upper>().solve(y2.topRows(rk));
186
187 // Apply the column permutation
188 // colsPermutation() performs a copy of the permutation,
189 // so let's apply it manually:
190 for(Index i = 0; i < rk; ++i) dest.row(m_E[i]) = y.row(i);
191 for(Index i = rk; i < cols(); ++i) dest.row(m_E[i]).setZero();
192
193// y.bottomRows(y.rows()-rk).setZero();
194// dest = colsPermutation() * y.topRows(cols());
195
196 m_info = Success;
197 }
198
201 const MatrixType matrixR() const
202 {
203 eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
204 if(!m_isRUpToDate) {
205 m_R = viewAsEigen<Scalar,ColMajor, typename MatrixType::StorageIndex>(*m_cR);
206 m_isRUpToDate = true;
207 }
208 return m_R;
209 }
217 {
218 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
219 return PermutationType(m_E, m_cR->ncol);
220 }
225 Index rank() const
226 {
227 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
228 return m_cc.SPQR_istat[4];
229 }
231 void setSPQROrdering(int ord) { m_ordering = ord;}
234 {
235 m_useDefaultThreshold = false;
236 m_tolerance = tol;
237 }
238
240 cholmod_common *cholmodCommon() const { return &m_cc; }
241
242
249 {
250 eigen_assert(m_isInitialized && "Decomposition is not initialized.");
251 return m_info;
252 }
253 protected:
256 mutable bool m_isRUpToDate;
258 int m_ordering; // Ordering method to use, see SPQR's manual
259 int m_allow_tol; // Allow to use some tolerance during numerical factorization.
260 RealScalar m_tolerance; // treat columns with 2-norm below this tolerance as zero
261 mutable cholmod_sparse *m_cR; // The sparse R factor in cholmod format
262 mutable MatrixType m_R; // The sparse matrix R in Eigen format
263 mutable StorageIndex *m_E; // The permutation applied to columns
264 mutable cholmod_sparse *m_H; //The householder vectors
265 mutable StorageIndex *m_HPinv; // The row permutation of H
266 mutable cholmod_dense *m_HTau; // The Householder coefficients
267 mutable Index m_rank; // The rank of the matrix
268 mutable cholmod_common m_cc; // Workspace and parameters
269 bool m_useDefaultThreshold; // Use default threshold
271 template<typename ,typename > friend struct SPQR_QProduct;
272};
273
274template <typename SPQRType, typename Derived>
275struct SPQR_QProduct : ReturnByValue<SPQR_QProduct<SPQRType,Derived> >
276{
277 typedef typename SPQRType::Scalar Scalar;
278 typedef typename SPQRType::StorageIndex StorageIndex;
279 //Define the constructor to get reference to argument types
280 SPQR_QProduct(const SPQRType& spqr, const Derived& other, bool transpose) : m_spqr(spqr),m_other(other),m_transpose(transpose) {}
281
282 inline Index rows() const { return m_transpose ? m_spqr.rows() : m_spqr.cols(); }
283 inline Index cols() const { return m_other.cols(); }
284 // Assign to a vector
285 template<typename ResType>
286 void evalTo(ResType& res) const
287 {
288 cholmod_dense y_cd;
289 cholmod_dense *x_cd;
290 int method = m_transpose ? SPQR_QTX : SPQR_QX;
291 cholmod_common *cc = m_spqr.cholmodCommon();
292 y_cd = viewAsCholmod(m_other.const_cast_derived());
293 x_cd = SuiteSparseQR_qmult<Scalar>(method, m_spqr.m_H, m_spqr.m_HTau, m_spqr.m_HPinv, &y_cd, cc);
294 res = Matrix<Scalar,ResType::RowsAtCompileTime,ResType::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x), x_cd->nrow, x_cd->ncol);
295 cholmod_l_free_dense(&x_cd, cc);
296 }
297 const SPQRType& m_spqr;
298 const Derived& m_other;
300
301};
302template<typename SPQRType>
304
305 SPQRMatrixQReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
306 template<typename Derived>
315 // To use for operations with the transpose of Q
320 const SPQRType& m_spqr;
321};
322
323template<typename SPQRType>
325 SPQRMatrixQTransposeReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
326 template<typename Derived>
328 {
329 return SPQR_QProduct<SPQRType,Derived>(m_spqr,other.derived(), true);
330 }
331 const SPQRType& m_spqr;
332};
333
334}// End namespace Eigen
335#endif
int i
Definition BiCGSTAB_step_by_step.cpp:9
#define eigen_assert(x)
Definition Macros.h:1037
m col(1)
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Definition PartialRedux_count.cpp:3
MatrixXf mat
Definition Tutorial_AdvancedInitialization_CommaTemporary.cpp:1
Scalar * b
Definition benchVecAdd.cpp:17
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition bench_gemm.cpp:48
EIGEN_DEVICE_FUNC Derived & setZero()
Definition CwiseNullaryOp.h:546
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:50
The matrix class, also used for vectors and row-vectors.
Definition Matrix.h:180
Definition ReturnByValue.h:52
Sparse QR factorization based on SuiteSparseQR library.
Definition SuiteSparseQRSupport.h:61
cholmod_sparse * m_H
Definition SuiteSparseQRSupport.h:264
StorageIndex * m_HPinv
Definition SuiteSparseQRSupport.h:265
MatrixType m_R
Definition SuiteSparseQRSupport.h:262
cholmod_common m_cc
Definition SuiteSparseQRSupport.h:268
bool m_analysisIsOk
Definition SuiteSparseQRSupport.h:254
ComputationInfo info() const
Reports whether previous computation was successful.
Definition SuiteSparseQRSupport.h:248
StorageIndex * m_E
Definition SuiteSparseQRSupport.h:263
Index m_rank
Definition SuiteSparseQRSupport.h:267
SPQR(const _MatrixType &matrix)
Definition SuiteSparseQRSupport.h:93
bool m_useDefaultThreshold
Definition SuiteSparseQRSupport.h:269
Index rank() const
Definition SuiteSparseQRSupport.h:225
cholmod_dense * m_HTau
Definition SuiteSparseQRSupport.h:266
SuiteSparse_long StorageIndex
Definition SuiteSparseQRSupport.h:68
SparseSolverBase< SPQR< _MatrixType > > Base
Definition SuiteSparseQRSupport.h:63
Index rows() const
Definition SuiteSparseQRSupport.h:164
void SPQR_free()
Definition SuiteSparseQRSupport.h:116
_MatrixType::Scalar Scalar
Definition SuiteSparseQRSupport.h:66
bool m_factorizationIsOk
Definition SuiteSparseQRSupport.h:255
ComputationInfo m_info
Definition SuiteSparseQRSupport.h:257
Index m_rows
Definition SuiteSparseQRSupport.h:270
int m_allow_tol
Definition SuiteSparseQRSupport.h:259
SPQR()
Definition SuiteSparseQRSupport.h:76
void compute(const _MatrixType &matrix)
Definition SuiteSparseQRSupport.h:125
Index cols() const
Definition SuiteSparseQRSupport.h:169
_MatrixType::RealScalar RealScalar
Definition SuiteSparseQRSupport.h:67
RealScalar m_tolerance
Definition SuiteSparseQRSupport.h:260
Map< PermutationMatrix< Dynamic, Dynamic, StorageIndex > > PermutationType
Definition SuiteSparseQRSupport.h:70
PermutationType colsPermutation() const
Get the permutation that was applied to columns of A.
Definition SuiteSparseQRSupport.h:216
SPQRMatrixQReturnType< SPQR > matrixQ() const
Get an expression of the matrix Q.
Definition SuiteSparseQRSupport.h:211
void setPivotThreshold(const RealScalar &tol)
Set the tolerance tol to treat columns with 2-norm < =tol as zero.
Definition SuiteSparseQRSupport.h:233
int m_ordering
Definition SuiteSparseQRSupport.h:258
SparseMatrix< Scalar, ColMajor, StorageIndex > MatrixType
Definition SuiteSparseQRSupport.h:69
cholmod_common * cholmodCommon() const
Definition SuiteSparseQRSupport.h:240
const MatrixType matrixR() const
Definition SuiteSparseQRSupport.h:201
void _solve_impl(const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
Definition SuiteSparseQRSupport.h:172
cholmod_sparse * m_cR
Definition SuiteSparseQRSupport.h:261
void setSPQROrdering(int ord)
Set the fill-reducing ordering method to be used.
Definition SuiteSparseQRSupport.h:231
~SPQR()
Definition SuiteSparseQRSupport.h:111
bool m_isRUpToDate
Definition SuiteSparseQRSupport.h:256
@ MaxColsAtCompileTime
Definition SuiteSparseQRSupport.h:73
@ ColsAtCompileTime
Definition SuiteSparseQRSupport.h:72
A versatible sparse matrix representation.
Definition SparseMatrix.h:98
A base class for sparse solvers.
Definition SparseSolverBase.h:68
bool m_isInitialized
Definition SparseSolverBase.h:119
Map< Matrix< T, Dynamic, Dynamic, ColMajor >, 0, OuterStride<> > matrix(T *data, int rows, int cols, int stride)
Definition common.h:110
ComputationInfo
Definition Constants.h:440
@ NumericalIssue
Definition Constants.h:444
@ Success
Definition Constants.h:442
Scalar * y
Definition level1_cplx_impl.h:124
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition MathFunctions.h:1091
Namespace containing all symbols from the Eigen library.
Definition bench_norm.cpp:85
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:74
cholmod_sparse viewAsCholmod(Ref< SparseMatrix< _Scalar, _Options, _StorageIndex > > mat)
Definition CholmodSupport.h:58
const int Dynamic
Definition Constants.h:22
Definition BandTriangularSolver.h:13
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:233
Definition SuiteSparseQRSupport.h:303
SPQRMatrixQReturnType(const SPQRType &spqr)
Definition SuiteSparseQRSupport.h:305
const SPQRType & m_spqr
Definition SuiteSparseQRSupport.h:320
SPQRMatrixQTransposeReturnType< SPQRType > adjoint() const
Definition SuiteSparseQRSupport.h:311
SPQRMatrixQTransposeReturnType< SPQRType > transpose() const
Definition SuiteSparseQRSupport.h:316
SPQR_QProduct< SPQRType, Derived > operator*(const MatrixBase< Derived > &other)
Definition SuiteSparseQRSupport.h:307
Definition SuiteSparseQRSupport.h:324
SPQRMatrixQTransposeReturnType(const SPQRType &spqr)
Definition SuiteSparseQRSupport.h:325
const SPQRType & m_spqr
Definition SuiteSparseQRSupport.h:331
SPQR_QProduct< SPQRType, Derived > operator*(const MatrixBase< Derived > &other)
Definition SuiteSparseQRSupport.h:327
Definition SuiteSparseQRSupport.h:276
bool m_transpose
Definition SuiteSparseQRSupport.h:299
SPQRType::Scalar Scalar
Definition SuiteSparseQRSupport.h:277
SPQR_QProduct(const SPQRType &spqr, const Derived &other, bool transpose)
Definition SuiteSparseQRSupport.h:280
void evalTo(ResType &res) const
Definition SuiteSparseQRSupport.h:286
const SPQRType & m_spqr
Definition SuiteSparseQRSupport.h:297
SPQRType::StorageIndex StorageIndex
Definition SuiteSparseQRSupport.h:278
Index cols() const
Definition SuiteSparseQRSupport.h:283
Index rows() const
Definition SuiteSparseQRSupport.h:282
const Derived & m_other
Definition SuiteSparseQRSupport.h:298
SPQRType::MatrixType ReturnType
Definition SuiteSparseQRSupport.h:23
SPQRType::MatrixType ReturnType
Definition SuiteSparseQRSupport.h:27
Derived::PlainObject ReturnType
Definition SuiteSparseQRSupport.h:31
Definition ForwardDeclarations.h:17
std::ptrdiff_t j
Definition tut_arithmetic_redux_minmax.cpp:2