aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseQR/SparseQR.h
diff options
context:
space:
mode:
authorGravatar Christoph Hertzberg <chtz@informatik.uni-bremen.de>2014-12-04 22:48:53 +0100
committerGravatar Christoph Hertzberg <chtz@informatik.uni-bremen.de>2014-12-04 22:48:53 +0100
commite8cdbedefb1913b5a0e2f2b7d38470f081cb8d29 (patch)
treeb64cb33df57f4cfcd87bf42643279629dc0900d3 /Eigen/src/SparseQR/SparseQR.h
parent6ccf97f3e6ce39c210e225ba7aae66da15b71660 (diff)
bug #877, bug #572: Introduce a global Index typedef. Rename Sparse*::Index to StorageIndex, make Dense*::StorageIndex an alias to DenseIndex. Overall this commit gets rid of all Index conversion warnings.
Diffstat (limited to 'Eigen/src/SparseQR/SparseQR.h')
-rw-r--r--Eigen/src/SparseQR/SparseQR.h46
1 files changed, 23 insertions, 23 deletions
diff --git a/Eigen/src/SparseQR/SparseQR.h b/Eigen/src/SparseQR/SparseQR.h
index 133211488..58bfc1cb4 100644
--- a/Eigen/src/SparseQR/SparseQR.h
+++ b/Eigen/src/SparseQR/SparseQR.h
@@ -21,7 +21,7 @@ namespace internal {
template <typename SparseQRType> struct traits<SparseQRMatrixQReturnType<SparseQRType> >
{
typedef typename SparseQRType::MatrixType ReturnType;
- typedef typename ReturnType::Index Index;
+ typedef typename ReturnType::StorageIndex StorageIndex;
typedef typename ReturnType::StorageKind StorageKind;
};
template <typename SparseQRType> struct traits<SparseQRMatrixQTransposeReturnType<SparseQRType> >
@@ -73,11 +73,11 @@ class SparseQR : public SparseSolverBase<SparseQR<_MatrixType,_OrderingType> >
typedef _OrderingType OrderingType;
typedef typename MatrixType::Scalar Scalar;
typedef typename MatrixType::RealScalar RealScalar;
- typedef typename MatrixType::Index Index;
- typedef SparseMatrix<Scalar,ColMajor,Index> QRMatrixType;
- typedef Matrix<Index, Dynamic, 1> IndexVector;
+ typedef typename MatrixType::StorageIndex StorageIndex;
+ typedef SparseMatrix<Scalar,ColMajor,StorageIndex> QRMatrixType;
+ typedef Matrix<StorageIndex, Dynamic, 1> IndexVector;
typedef Matrix<Scalar, Dynamic, 1> ScalarVector;
- typedef PermutationMatrix<Dynamic, Dynamic, Index> PermutationType;
+ typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> PermutationType;
public:
SparseQR () : m_analysisIsok(false), m_lastError(""), m_useDefaultThreshold(true),m_isQSorted(false),m_isEtreeOk(false)
{ }
@@ -109,11 +109,11 @@ class SparseQR : public SparseSolverBase<SparseQR<_MatrixType,_OrderingType> >
/** \returns the number of rows of the represented matrix.
*/
- inline Index rows() const { return m_pmat.rows(); }
+ inline StorageIndex rows() const { return m_pmat.rows(); }
/** \returns the number of columns of the represented matrix.
*/
- inline Index cols() const { return m_pmat.cols();}
+ inline StorageIndex cols() const { return m_pmat.cols();}
/** \returns a const reference to the \b sparse upper triangular matrix R of the QR factorization.
*/
@@ -123,7 +123,7 @@ class SparseQR : public SparseSolverBase<SparseQR<_MatrixType,_OrderingType> >
*
* \sa setPivotThreshold()
*/
- Index rank() const
+ StorageIndex rank() const
{
eigen_assert(m_isInitialized && "The factorization should be called first, use compute()");
return m_nonzeropivots;
@@ -179,7 +179,7 @@ class SparseQR : public SparseSolverBase<SparseQR<_MatrixType,_OrderingType> >
b = y;
// Solve with the triangular matrix R
- y.resize((std::max)(cols(),Index(y.rows())),y.cols());
+ y.resize((std::max<Index>)(cols(),y.rows()),y.cols());
y.topRows(rank) = this->matrixR().topLeftCorner(rank, rank).template triangularView<Upper>().solve(b.topRows(rank));
y.bottomRows(y.rows()-rank).setZero();
@@ -260,7 +260,7 @@ class SparseQR : public SparseSolverBase<SparseQR<_MatrixType,_OrderingType> >
PermutationType m_outputPerm_c; // The final column permutation
RealScalar m_threshold; // Threshold to determine null Householder reflections
bool m_useDefaultThreshold; // Use default threshold
- Index m_nonzeropivots; // Number of non zero pivots found
+ StorageIndex m_nonzeropivots; // Number of non zero pivots found
IndexVector m_etree; // Column elimination tree
IndexVector m_firstRowElt; // First element in each row
bool m_isQSorted; // whether Q is sorted or not
@@ -289,9 +289,9 @@ void SparseQR<MatrixType,OrderingType>::analyzePattern(const MatrixType& mat)
// Compute the column fill reducing ordering
OrderingType ord;
ord(matCpy, m_perm_c);
- Index n = mat.cols();
- Index m = mat.rows();
- Index diagSize = (std::min)(m,n);
+ StorageIndex n = mat.cols();
+ StorageIndex m = mat.rows();
+ StorageIndex diagSize = (std::min)(m,n);
if (!m_perm_c.size())
{
@@ -354,7 +354,7 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
// otherwise directly use the input matrix
//
IndexVector originalOuterIndicesCpy;
- const Index *originalOuterIndices = mat.outerIndexPtr();
+ const StorageIndex *originalOuterIndices = mat.outerIndexPtr();
if(MatrixType::IsRowMajor)
{
originalOuterIndicesCpy = IndexVector::Map(m_pmat.outerIndexPtr(),n+1);
@@ -385,11 +385,11 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
// Initialize the numerical permutation
m_pivotperm.setIdentity(n);
- Index nonzeroCol = 0; // Record the number of valid pivots
+ StorageIndex nonzeroCol = 0; // Record the number of valid pivots
m_Q.startVec(0);
// Left looking rank-revealing QR factorization: compute a column of R and Q at a time
- for (Index col = 0; col < n; ++col)
+ for (StorageIndex col = 0; col < n; ++col)
{
mark.setConstant(-1);
m_R.startVec(col);
@@ -405,12 +405,12 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
// thus the trick with found_diag that permits to do one more iteration on the diagonal element if this one has not been found.
for (typename QRMatrixType::InnerIterator itp(m_pmat, col); itp || !found_diag; ++itp)
{
- Index curIdx = nonzeroCol;
+ StorageIndex curIdx = nonzeroCol;
if(itp) curIdx = itp.row();
if(curIdx == nonzeroCol) found_diag = true;
// Get the nonzeros indexes of the current column of R
- Index st = m_firstRowElt(curIdx); // The traversal of the etree starts here
+ StorageIndex st = m_firstRowElt(curIdx); // The traversal of the etree starts here
if (st < 0 )
{
m_lastError = "Empty row found during numerical factorization";
@@ -467,7 +467,7 @@ void SparseQR<MatrixType,OrderingType>::factorize(const MatrixType& mat)
{
for (typename QRMatrixType::InnerIterator itq(m_Q, curIdx); itq; ++itq)
{
- Index iQ = itq.row();
+ StorageIndex iQ = itq.row();
if (mark(iQ) != col)
{
Qidx(nzcolQ++) = iQ; // Add this row to the pattern of Q,
@@ -578,7 +578,7 @@ struct SparseQR_QProduct : ReturnByValue<SparseQR_QProduct<SparseQRType, Derived
{
typedef typename SparseQRType::QRMatrixType MatrixType;
typedef typename SparseQRType::Scalar Scalar;
- typedef typename SparseQRType::Index Index;
+ typedef typename SparseQRType::StorageIndex StorageIndex;
// Get the references
SparseQR_QProduct(const SparseQRType& qr, const Derived& other, bool transpose) :
m_qr(qr),m_other(other),m_transpose(transpose) {}
@@ -634,7 +634,7 @@ struct SparseQR_QProduct : ReturnByValue<SparseQR_QProduct<SparseQRType, Derived
template<typename SparseQRType>
struct SparseQRMatrixQReturnType : public EigenBase<SparseQRMatrixQReturnType<SparseQRType> >
{
- typedef typename SparseQRType::Index Index;
+ typedef typename SparseQRType::StorageIndex StorageIndex;
typedef typename SparseQRType::Scalar Scalar;
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
explicit SparseQRMatrixQReturnType(const SparseQRType& qr) : m_qr(qr) {}
@@ -647,8 +647,8 @@ struct SparseQRMatrixQReturnType : public EigenBase<SparseQRMatrixQReturnType<Sp
{
return SparseQRMatrixQTransposeReturnType<SparseQRType>(m_qr);
}
- inline Index rows() const { return m_qr.rows(); }
- inline Index cols() const { return (std::min)(m_qr.rows(),m_qr.cols()); }
+ inline StorageIndex rows() const { return m_qr.rows(); }
+ inline StorageIndex cols() const { return (std::min)(m_qr.rows(),m_qr.cols()); }
// To use for operations with the transpose of Q
SparseQRMatrixQTransposeReturnType<SparseQRType> transpose() const
{