11 #ifndef EIGEN_SPARSE_QR_H 12 #define EIGEN_SPARSE_QR_H 16 template<
typename MatrixType,
typename OrderingType>
class SparseQR;
71 template<
typename _MatrixType,
typename _OrderingType>
76 using Base::m_isInitialized;
78 using Base::_solve_impl;
90 ColsAtCompileTime = MatrixType::ColsAtCompileTime,
91 MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
95 SparseQR () : m_analysisIsok(false), m_lastError(
""), m_useDefaultThreshold(true),m_isQSorted(false),m_isEtreeOk(false)
104 explicit SparseQR(
const MatrixType&
mat) : m_analysisIsok(false), m_lastError(
""), m_useDefaultThreshold(true),m_isQSorted(false),m_isEtreeOk(false)
120 void analyzePattern(
const MatrixType&
mat);
121 void factorize(
const MatrixType& mat);
144 const QRMatrixType&
matrixR()
const {
return m_R; }
152 eigen_assert(m_isInitialized &&
"The factorization should be called first, use compute()");
153 return m_nonzeropivots;
182 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
183 return m_outputPerm_c;
192 template<
typename Rhs,
typename Dest>
195 eigen_assert(m_isInitialized &&
"The factorization should be called first, use compute()");
196 eigen_assert(this->
rows() == B.rows() &&
"SparseQR::solve() : invalid number of rows in the right hand side matrix");
198 Index rank = this->rank();
201 typename Dest::PlainObject
y,
b;
202 y = this->matrixQ().adjoint() * B;
206 y.resize((std::max<Index>)(
cols(),y.rows()),y.cols());
207 y.topRows(rank) = this->matrixR().topLeftCorner(rank, rank).template triangularView<Upper>().solve(b.topRows(rank));
208 y.bottomRows(y.rows()-rank).
setZero();
211 if (m_perm_c.size()) dest = colsPermutation() * y.
topRows(
cols());
225 m_useDefaultThreshold =
false;
226 m_threshold = threshold;
233 template<
typename Rhs>
236 eigen_assert(m_isInitialized &&
"The factorization should be called first, use compute()");
237 eigen_assert(this->
rows() == B.rows() &&
"SparseQR::solve() : invalid number of rows in the right hand side matrix");
240 template<
typename Rhs>
243 eigen_assert(m_isInitialized &&
"The factorization should be called first, use compute()");
244 eigen_assert(this->
rows() == B.
rows() &&
"SparseQR::solve() : invalid number of rows in the right hand side matrix");
258 eigen_assert(m_isInitialized &&
"Decomposition is not initialized.");
266 if(this->m_isQSorted)
return;
270 this->m_isQSorted =
true;
307 template <
typename MatrixType,
typename OrderingType>
310 eigen_assert(mat.isCompressed() &&
"SparseQR requires a sparse matrix in compressed mode. Call .makeCompressed() before passing it to SparseQR");
315 ord(matCpy, m_perm_c);
320 if (!m_perm_c.size())
323 m_perm_c.indices().setLinSpaced(n, 0,
StorageIndex(n-1));
327 m_outputPerm_c = m_perm_c.inverse();
332 m_Q.resize(m, diagSize);
335 m_R.reserve(2*mat.nonZeros());
336 m_Q.reserve(2*mat.nonZeros());
337 m_hcoeffs.resize(diagSize);
338 m_analysisIsok =
true;
348 template <
typename MatrixType,
typename OrderingType>
353 eigen_assert(m_analysisIsok &&
"analyzePattern() should be called before this step");
359 Index nzcolR, nzcolQ;
368 m_outputPerm_c = m_perm_c.inverse();
381 const StorageIndex *originalOuterIndices = mat.outerIndexPtr();
382 if(MatrixType::IsRowMajor)
384 originalOuterIndicesCpy = IndexVector::Map(m_pmat.outerIndexPtr(),n+1);
385 originalOuterIndices = originalOuterIndicesCpy.
data();
388 for (
int i = 0;
i <
n;
i++)
390 Index p = m_perm_c.size() ? m_perm_c.indices()(
i) :
i;
391 m_pmat.outerIndexPtr()[
p] = originalOuterIndices[
i];
392 m_pmat.innerNonZeroPtr()[
p] = originalOuterIndices[
i+1] - originalOuterIndices[
i];
400 if(m_useDefaultThreshold)
403 for (
int j = 0;
j <
n;
j++) max2Norm =
numext::maxi(max2Norm, m_pmat.col(
j).norm());
410 m_pivotperm.setIdentity(n);
420 mark(nonzeroCol) =
col;
421 Qidx(0) = nonzeroCol;
422 nzcolR = 0; nzcolQ = 1;
423 bool found_diag = nonzeroCol>=
m;
434 if(curIdx == nonzeroCol) found_diag =
true;
440 m_lastError =
"Empty row found during numerical factorization";
447 for (; mark(st) !=
col; st = m_etree(st))
455 Index nt = nzcolR-bi;
459 if(itp) tval(curIdx) = itp.value();
460 else tval(curIdx) =
Scalar(0);
463 if(curIdx > nonzeroCol && mark(curIdx) !=
col )
465 Qidx(nzcolQ) = curIdx;
472 for (
Index i = nzcolR-1;
i >= 0;
i--)
480 tdot = m_Q.col(curIdx).dot(tval);
482 tdot *= m_hcoeffs(curIdx);
487 tval(itq.row()) -= itq.value() * tdot;
490 if(m_etree(Ridx(
i)) == nonzeroCol)
507 if(nonzeroCol < diagSize)
515 for (
Index itq = 1; itq < nzcolQ; ++itq) sqrNorm +=
numext::abs2(tval(Qidx(itq)));
528 for (
Index itq = 1; itq < nzcolQ; ++itq)
529 tval(Qidx(itq)) /= (c0 - beta);
536 for (
Index i = nzcolR-1;
i >= 0;
i--)
539 if(curIdx < nonzeroCol)
541 m_R.insertBackByOuterInnerUnordered(
col, curIdx) = tval(curIdx);
542 tval(curIdx) =
Scalar(0.);
546 if(nonzeroCol < diagSize &&
abs(beta) >= pivotThreshold)
548 m_R.insertBackByOuterInner(
col, nonzeroCol) = beta;
550 m_hcoeffs(nonzeroCol) = tau;
552 for (
Index itq = 0; itq < nzcolQ; ++itq)
554 Index iQ = Qidx(itq);
555 m_Q.insertBackByOuterInnerUnordered(nonzeroCol,iQ) = tval(iQ);
559 if(nonzeroCol<diagSize)
560 m_Q.startVec(nonzeroCol);
565 for (
Index j = nonzeroCol;
j < n-1;
j++)
566 std::swap(m_pivotperm.indices()(
j), m_pivotperm.indices()[
j+1]);
574 m_hcoeffs.tail(diagSize-nonzeroCol).setZero();
578 m_Q.makeCompressed();
580 m_R.makeCompressed();
583 m_nonzeropivots = nonzeroCol;
589 m_R = tempR * m_pivotperm;
592 m_outputPerm_c = m_outputPerm_c * m_pivotperm;
595 m_isInitialized =
true;
596 m_factorizationIsok =
true;
600 template <
typename SparseQRType,
typename Derived>
607 m_qr(qr),m_other(other),m_transpose(transpose) {}
608 inline Index rows()
const {
return m_qr.matrixQ().rows(); }
612 template<
typename DesType>
621 eigen_assert(m_qr.m_Q.rows() == m_other.rows() &&
"Non conforming object sizes");
623 for(
Index j = 0;
j < res.cols();
j++){
624 for (
Index k = 0; k < diagSize; k++)
627 tau = m_qr.m_Q.col(k).dot(res.col(
j));
628 if(tau==
Scalar(0))
continue;
629 tau = tau * m_qr.m_hcoeffs(k);
630 res.col(
j) -= tau * m_qr.m_Q.col(k);
636 eigen_assert(m_qr.matrixQ().cols() == m_other.rows() &&
"Non conforming object sizes");
638 res.conservativeResize(
rows(),
cols());
641 for(
Index j = 0;
j < res.cols();
j++)
643 for (
Index k = diagSize-1; k >=0; k--)
646 tau = m_qr.m_Q.col(k).dot(res.col(
j));
647 if(tau==
Scalar(0))
continue;
649 res.col(
j) -= tau * m_qr.m_Q.col(k);
660 template<
typename SparseQRType>
670 template<
typename Derived>
691 template<
typename SparseQRType>
695 template<
typename Derived>
705 template<
typename SparseQRType>
713 template<
typename DstXprType,
typename SparseQRType>
721 typename DstXprType::PlainObject idMat(src.
rows(), src.
cols());
724 const_cast<SparseQRType *
>(&src.
m_qr)->_sort_matrix_Q();
729 template<
typename DstXprType,
typename SparseQRType>
737 dst = src.
m_qr.matrixQ() * DstXprType::Identity(src.
m_qr.rows(), src.
m_qr.rows());
_OrderingType OrderingType
MatrixType::StorageIndex StorageIndex
PermutationMatrix< Dynamic, Dynamic, StorageIndex > PermutationType
Matrix< Scalar, Dynamic, 1 > ScalarVector
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op< Scalar, Scalar > &)
EIGEN_DEVICE_FUNC Derived & setZero(Index size)
Derived::PlainObject ReturnType
int coletree(const MatrixType &mat, IndexVector &parent, IndexVector &firstRowElt, typename MatrixType::StorageIndex *perm=0)
const SparseQRType & m_qr
SparseQRType::MatrixType ReturnType
SparseQRMatrixQReturnType(const SparseQRType &qr)
Matrix< StorageIndex, Dynamic, 1 > IndexVector
const SparseQRType & m_qr
A base class for sparse solvers.
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Scalar * data() const
storage_kind_to_evaluator_kind< typename MatrixType::StorageKind >::Kind Kind
SparseQR_QProduct(const SparseQRType &qr, const Derived &other, bool transpose)
EIGEN_DEVICE_FUNC const SqrtReturnType sqrt() const
Namespace containing all symbols from the Eigen library.
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Eigen::Index Index
The interface type of indices.
SparseQRMatrixQReturnType< SparseQRType > SrcXprType
SparseQRType::Scalar Scalar
EIGEN_DEVICE_FUNC RowsBlockXpr topRows(Index n)
ComputationInfo info() const
Reports whether previous computation was successful.
const QRMatrixType & matrixR() const
DstXprType::Scalar Scalar
SparseQRType::QRMatrixType MatrixType
HouseholderQR< MatrixXf > qr(A)
void factorize(const MatrixType &mat)
Performs the numerical QR factorization of the input matrix.
SparseQR_QProduct< SparseQRType, Derived > operator*(const MatrixBase< Derived > &other)
static void run(DstXprType &dst, const SrcXprType &src, const internal::assign_op< Scalar, Scalar > &)
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
SparseQRType::MatrixType MatrixType
Sparse left-looking rank-revealing QR factorization.
IndexVector m_firstRowElt
Base class of any sparse matrices or sparse expressions.
SparseSolverBase< SparseQR< _MatrixType, _OrderingType > > Base
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
Base::InnerIterator InnerIterator
PermutationType m_pivotperm
DstXprType::Scalar Scalar
SparseQRMatrixQTransposeReturnType< SparseQRType > adjoint() const
EIGEN_DEVICE_FUNC Derived & setConstant(Index size, const Scalar &val)
const Solve< SparseQR, Rhs > solve(const SparseMatrixBase< Rhs > &B) const
NumTraits< Scalar >::Real RealScalar
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE const Abs2ReturnType abs2() const
MatrixType::Scalar Scalar
void evalTo(DesType &res) const
bool _solve_impl(const MatrixBase< Rhs > &B, MatrixBase< Dest > &dest) const
void analyzePattern(const MatrixType &mat)
Preprocessing step of a QR factorization.
SparseQRType::Scalar Scalar
PermutationType m_outputPerm_c
SparseQRMatrixQReturnType< SparseQRType > SrcXprType
SparseQRType::MatrixType ReturnType
const Derived & derived() const
const PermutationType & colsPermutation() const
const Solve< SparseQR, Rhs > solve(const MatrixBase< Rhs > &B) const
DstXprType::StorageIndex StorageIndex
SparseQRMatrixQReturnType< SparseQR > matrixQ() const
const SparseQRType & m_qr
EIGEN_DEVICE_FUNC const ImagReturnType imag() const
SparseQRMatrixQTransposeReturnType(const SparseQRType &qr)
ReturnType::StorageIndex StorageIndex
Pseudo expression representing a solving operation.
MatrixType::RealScalar RealScalar
EIGEN_DONT_INLINE void compute(Solver &solver, const MatrixType &A)
ReturnType::StorageKind StorageKind
SparseQR(const MatrixType &mat)
Base class for all dense matrices, vectors, and expressions.
void swap(mpfr::mpreal &x, mpfr::mpreal &y)
SparseQRMatrixQTransposeReturnType< SparseQRType > transpose() const
DstXprType::StorageIndex StorageIndex
void setPivotThreshold(const RealScalar &threshold)
std::string lastErrorMessage() const
bool m_useDefaultThreshold
ScalarWithExceptions conj(const ScalarWithExceptions &x)
SparseMatrix< Scalar, ColMajor, StorageIndex > QRMatrixType
SparseQR_QProduct< SparseQRType, Derived > operator*(const MatrixBase< Derived > &other)
void compute(const MatrixType &mat)