From fc202bab398ed9b78ef8452f8e4ef35e233965f6 Mon Sep 17 00:00:00 2001 From: Gael Guennebaud Date: Fri, 13 Feb 2015 18:57:41 +0100 Subject: Index refactoring: StorageIndex must be used for storage only (and locally when it make sense). In all other cases use the global Index type. --- Eigen/src/SparseLU/SparseLU.h | 26 ++++++++++++-------------- Eigen/src/SparseLU/SparseLUImpl.h | 8 ++++---- Eigen/src/SparseLU/SparseLU_Memory.h | 15 +++++++-------- Eigen/src/SparseLU/SparseLU_Structs.h | 3 +-- Eigen/src/SparseLU/SparseLU_SupernodalMatrix.h | 18 +++++++++--------- Eigen/src/SparseLU/SparseLU_Utils.h | 8 ++++---- Eigen/src/SparseLU/SparseLU_column_bmod.h | 5 +++-- Eigen/src/SparseLU/SparseLU_column_dfs.h | 16 +++++++++------- Eigen/src/SparseLU/SparseLU_copy_to_ucol.h | 5 +++-- Eigen/src/SparseLU/SparseLU_gemm_kernel.h | 2 +- Eigen/src/SparseLU/SparseLU_heap_relax_snode.h | 6 +++--- Eigen/src/SparseLU/SparseLU_kernel_bmod.h | 11 ++++++----- Eigen/src/SparseLU/SparseLU_panel_bmod.h | 4 ++-- Eigen/src/SparseLU/SparseLU_panel_dfs.h | 14 +++++++------- Eigen/src/SparseLU/SparseLU_pivotL.h | 6 +++--- Eigen/src/SparseLU/SparseLU_pruneL.h | 5 +++-- Eigen/src/SparseLU/SparseLU_relax_snode.h | 4 ++-- 17 files changed, 79 insertions(+), 77 deletions(-) (limited to 'Eigen/src/SparseLU') diff --git a/Eigen/src/SparseLU/SparseLU.h b/Eigen/src/SparseLU/SparseLU.h index 380ba25c0..f60e4ac9d 100644 --- a/Eigen/src/SparseLU/SparseLU.h +++ b/Eigen/src/SparseLU/SparseLU.h @@ -122,8 +122,8 @@ class SparseLU : public SparseSolverBase >, factorize(matrix); } - inline StorageIndex rows() const { return m_mat.rows(); } - inline StorageIndex cols() const { return m_mat.cols(); } + inline Index rows() const { return m_mat.rows(); } + inline Index cols() const { return m_mat.cols(); } /** Indicate that the pattern of the input matrix is symmetric */ void isSymmetric(bool sym) { @@ -334,10 +334,10 @@ class SparseLU : public SparseSolverBase >, // SparseLU options bool m_symmetricmode; // values for performance - internal::perfvalues m_perfv; + internal::perfvalues m_perfv; RealScalar m_diagpivotthresh; // Specifies the threshold used for a diagonal entry to be an acceptable pivot - StorageIndex m_nnzL, m_nnzU; // Nonzeros in L and U factors - StorageIndex m_detPermR; // Determinant of the coefficient matrix + Index m_nnzL, m_nnzU; // Nonzeros in L and U factors + Index m_detPermR; // Determinant of the coefficient matrix private: // Disable copy constructor SparseLU (const SparseLU& ); @@ -449,7 +449,7 @@ void SparseLU::factorize(const MatrixType& matrix) eigen_assert(m_analysisIsOk && "analyzePattern() should be called first"); eigen_assert((matrix.rows() == matrix.cols()) && "Only for squared matrices"); - typedef typename IndexVector::Scalar Index; + typedef typename IndexVector::Scalar StorageIndex; m_isInitialized = true; @@ -461,11 +461,11 @@ void SparseLU::factorize(const MatrixType& matrix) { m_mat.uncompress(); //NOTE: The effect of this command is only to create the InnerNonzeros pointers. //Then, permute only the column pointers - const Index * outerIndexPtr; + const StorageIndex * outerIndexPtr; if (matrix.isCompressed()) outerIndexPtr = matrix.outerIndexPtr(); else { - Index* outerIndexPtr_t = new Index[matrix.cols()+1]; + StorageIndex* outerIndexPtr_t = new StorageIndex[matrix.cols()+1]; for(Index i = 0; i <= matrix.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i]; outerIndexPtr = outerIndexPtr_t; } @@ -649,12 +649,11 @@ void SparseLU::factorize(const MatrixType& matrix) template struct SparseLUMatrixLReturnType : internal::no_assignment_operator { - typedef typename MappedSupernodalType::StorageIndex StorageIndex; typedef typename MappedSupernodalType::Scalar Scalar; explicit SparseLUMatrixLReturnType(const MappedSupernodalType& mapL) : m_mapL(mapL) { } - StorageIndex rows() { return m_mapL.rows(); } - StorageIndex cols() { return m_mapL.cols(); } + Index rows() { return m_mapL.rows(); } + Index cols() { return m_mapL.cols(); } template void solveInPlace( MatrixBase &X) const { @@ -666,13 +665,12 @@ struct SparseLUMatrixLReturnType : internal::no_assignment_operator template struct SparseLUMatrixUReturnType : internal::no_assignment_operator { - typedef typename MatrixLType::StorageIndex StorageIndex; typedef typename MatrixLType::Scalar Scalar; explicit SparseLUMatrixUReturnType(const MatrixLType& mapL, const MatrixUType& mapU) : m_mapL(mapL),m_mapU(mapU) { } - StorageIndex rows() { return m_mapL.rows(); } - StorageIndex cols() { return m_mapL.cols(); } + Index rows() { return m_mapL.rows(); } + Index cols() { return m_mapL.cols(); } template void solveInPlace(MatrixBase &X) const { diff --git a/Eigen/src/SparseLU/SparseLUImpl.h b/Eigen/src/SparseLU/SparseLUImpl.h index 14d70897d..e735fd5c8 100644 --- a/Eigen/src/SparseLU/SparseLUImpl.h +++ b/Eigen/src/SparseLU/SparseLUImpl.h @@ -16,17 +16,17 @@ namespace internal { * \class SparseLUImpl * Base class for sparseLU */ -template +template class SparseLUImpl { public: typedef Matrix ScalarVector; - typedef Matrix IndexVector; + typedef Matrix IndexVector; typedef typename ScalarVector::RealScalar RealScalar; typedef Ref > BlockScalarVector; - typedef Ref > BlockIndexVector; + typedef Ref > BlockIndexVector; typedef LU_GlobalLU_t GlobalLU_t; - typedef SparseMatrix MatrixType; + typedef SparseMatrix MatrixType; protected: template diff --git a/Eigen/src/SparseLU/SparseLU_Memory.h b/Eigen/src/SparseLU/SparseLU_Memory.h index 1ffa7d54e..1cf8bebc7 100644 --- a/Eigen/src/SparseLU/SparseLU_Memory.h +++ b/Eigen/src/SparseLU/SparseLU_Memory.h @@ -36,13 +36,12 @@ namespace internal { enum { LUNoMarker = 3 }; enum {emptyIdxLU = -1}; -template inline Index LUnumTempV(Index& m, Index& w, Index& t, Index& b) { return (std::max)(m, (t+b)*w); } -template< typename Scalar, typename Index> +template< typename Scalar> inline Index LUTempSpace(Index&m, Index& w) { return (2*w + 4 + LUNoMarker) * m * sizeof(Index) + (w + 1) * m * sizeof(Scalar); @@ -59,9 +58,9 @@ inline Index LUTempSpace(Index&m, Index& w) * \param keep_prev 1: use length and do not expand the vector; 0: compute new_len and expand * \param[in,out] num_expansions Number of times the memory has been expanded */ -template +template template -Index SparseLUImpl::expand(VectorType& vec, Index& length, Index nbElts, Index keep_prev, Index& num_expansions) +Index SparseLUImpl::expand(VectorType& vec, Index& length, Index nbElts, Index keep_prev, Index& num_expansions) { float alpha = 1.5; // Ratio of the memory increase @@ -148,8 +147,8 @@ Index SparseLUImpl::expand(VectorType& vec, Index& length, Index * \return an estimated size of the required memory if lwork = -1; otherwise, return the size of actually allocated memory when allocation failed, and 0 on success * \note Unlike SuperLU, this routine does not support successive factorization with the same pattern and the same row permutation */ -template -Index SparseLUImpl::memInit(Index m, Index n, Index annz, Index lwork, Index fillratio, Index panel_size, GlobalLU_t& glu) +template +Index SparseLUImpl::memInit(Index m, Index n, Index annz, Index lwork, Index fillratio, Index panel_size, GlobalLU_t& glu) { Index& num_expansions = glu.num_expansions; //No memory expansions so far num_expansions = 0; @@ -205,9 +204,9 @@ Index SparseLUImpl::memInit(Index m, Index n, Index annz, Index lw * \param num_expansions Number of expansions * \return 0 on success, > 0 size of the memory allocated so far */ -template +template template -Index SparseLUImpl::memXpand(VectorType& vec, Index& maxlen, Index nbElts, MemType memtype, Index& num_expansions) +Index SparseLUImpl::memXpand(VectorType& vec, Index& maxlen, Index nbElts, MemType memtype, Index& num_expansions) { Index failed_size; if (memtype == USUB) diff --git a/Eigen/src/SparseLU/SparseLU_Structs.h b/Eigen/src/SparseLU/SparseLU_Structs.h index 24d6bf179..cf5ec449b 100644 --- a/Eigen/src/SparseLU/SparseLU_Structs.h +++ b/Eigen/src/SparseLU/SparseLU_Structs.h @@ -75,7 +75,7 @@ typedef enum {LUSUP, UCOL, LSUB, USUB, LLVL, ULVL} MemType; template struct LU_GlobalLU_t { - typedef typename IndexVector::Scalar Index; + typedef typename IndexVector::Scalar StorageIndex; IndexVector xsup; //First supernode column ... xsup(s) points to the beginning of the s-th supernode IndexVector supno; // Supernode number corresponding to this column (column to supernode mapping) ScalarVector lusup; // nonzero values of L ordered by columns @@ -93,7 +93,6 @@ struct LU_GlobalLU_t { }; // Values to set for performance -template struct perfvalues { Index panel_size; // a panel consists of at most consecutive columns Index relax; // To control degree of relaxing supernodes. If the number of nodes (columns) diff --git a/Eigen/src/SparseLU/SparseLU_SupernodalMatrix.h b/Eigen/src/SparseLU/SparseLU_SupernodalMatrix.h index 098763765..f7ffc2d9c 100644 --- a/Eigen/src/SparseLU/SparseLU_SupernodalMatrix.h +++ b/Eigen/src/SparseLU/SparseLU_SupernodalMatrix.h @@ -42,7 +42,7 @@ class MappedSuperNodalMatrix { } - MappedSuperNodalMatrix(StorageIndex m, StorageIndex n, ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind, + MappedSuperNodalMatrix(Index m, Index n, ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind, IndexVector& rowind_colptr, IndexVector& col_to_sup, IndexVector& sup_to_col ) { setInfos(m, n, nzval, nzval_colptr, rowind, rowind_colptr, col_to_sup, sup_to_col); @@ -58,7 +58,7 @@ class MappedSuperNodalMatrix * FIXME This class will be modified such that it can be use in the course * of the factorization. */ - void setInfos(StorageIndex m, StorageIndex n, ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind, + void setInfos(Index m, Index n, ScalarVector& nzval, IndexVector& nzval_colptr, IndexVector& rowind, IndexVector& rowind_colptr, IndexVector& col_to_sup, IndexVector& sup_to_col ) { m_row = m; @@ -75,12 +75,12 @@ class MappedSuperNodalMatrix /** * Number of rows */ - StorageIndex rows() { return m_row; } + Index rows() { return m_row; } /** * Number of columns */ - StorageIndex cols() { return m_col; } + Index cols() { return m_col; } /** * Return the array of nonzero values packed by column @@ -148,7 +148,7 @@ class MappedSuperNodalMatrix /** * Return the number of supernodes */ - StorageIndex nsuper() const + Index nsuper() const { return m_nsuper; } @@ -161,9 +161,9 @@ class MappedSuperNodalMatrix protected: - StorageIndex m_row; // Number of rows - StorageIndex m_col; // Number of columns - StorageIndex m_nsuper; // Number of supernodes + Index m_row; // Number of rows + Index m_col; // Number of columns + Index m_nsuper; // Number of supernodes Scalar* m_nzval; //array of nonzero values packed by column StorageIndex* m_nzval_colptr; //nzval_colptr[j] Stores the location in nzval[] which starts column j StorageIndex* m_rowind; // Array of compressed row indices of rectangular supernodes @@ -184,7 +184,7 @@ class MappedSuperNodalMatrix::InnerIterator public: InnerIterator(const MappedSuperNodalMatrix& mat, Eigen::Index outer) : m_matrix(mat), - m_outer(convert_index(outer)), + m_outer(outer), m_supno(mat.colToSup()[outer]), m_idval(mat.colIndexPtr()[outer]), m_startidval(m_idval), diff --git a/Eigen/src/SparseLU/SparseLU_Utils.h b/Eigen/src/SparseLU/SparseLU_Utils.h index 15352ac33..b48157d9f 100644 --- a/Eigen/src/SparseLU/SparseLU_Utils.h +++ b/Eigen/src/SparseLU/SparseLU_Utils.h @@ -17,8 +17,8 @@ namespace internal { /** * \brief Count Nonzero elements in the factors */ -template -void SparseLUImpl::countnz(const Index n, Index& nnzL, Index& nnzU, GlobalLU_t& glu) +template +void SparseLUImpl::countnz(const Index n, Index& nnzL, Index& nnzU, GlobalLU_t& glu) { nnzL = 0; nnzU = (glu.xusub)(n); @@ -48,8 +48,8 @@ void SparseLUImpl::countnz(const Index n, Index& nnzL, Index& nnzU * and applies permutation to the remaining subscripts * */ -template -void SparseLUImpl::fixupL(const Index n, const IndexVector& perm_r, GlobalLU_t& glu) +template +void SparseLUImpl::fixupL(const Index n, const IndexVector& perm_r, GlobalLU_t& glu) { Index fsupc, i, j, k, jstart; diff --git a/Eigen/src/SparseLU/SparseLU_column_bmod.h b/Eigen/src/SparseLU/SparseLU_column_bmod.h index f24bd87d3..bda01dcb3 100644 --- a/Eigen/src/SparseLU/SparseLU_column_bmod.h +++ b/Eigen/src/SparseLU/SparseLU_column_bmod.h @@ -49,8 +49,9 @@ namespace internal { * > 0 - number of bytes allocated when run out of space * */ -template -Index SparseLUImpl::column_bmod(const Index jcol, const Index nseg, BlockScalarVector dense, ScalarVector& tempv, BlockIndexVector segrep, BlockIndexVector repfnz, Index fpanelc, GlobalLU_t& glu) +template +Index SparseLUImpl::column_bmod(const Index jcol, const Index nseg, BlockScalarVector dense, ScalarVector& tempv, + BlockIndexVector segrep, BlockIndexVector repfnz, Index fpanelc, GlobalLU_t& glu) { Index jsupno, k, ksub, krep, ksupno; Index lptr, nrow, isub, irow, nextlu, new_next, ufirst; diff --git a/Eigen/src/SparseLU/SparseLU_column_dfs.h b/Eigen/src/SparseLU/SparseLU_column_dfs.h index 4c04b0e44..17c9e6adb 100644 --- a/Eigen/src/SparseLU/SparseLU_column_dfs.h +++ b/Eigen/src/SparseLU/SparseLU_column_dfs.h @@ -30,7 +30,7 @@ #ifndef SPARSELU_COLUMN_DFS_H #define SPARSELU_COLUMN_DFS_H -template class SparseLUImpl; +template class SparseLUImpl; namespace Eigen { namespace internal { @@ -39,8 +39,8 @@ template struct column_dfs_traits : no_assignment_operator { typedef typename ScalarVector::Scalar Scalar; - typedef typename IndexVector::Scalar Index; - column_dfs_traits(Index jcol, Index& jsuper, typename SparseLUImpl::GlobalLU_t& glu, SparseLUImpl& luImpl) + typedef typename IndexVector::Scalar StorageIndex; + column_dfs_traits(Index jcol, Index& jsuper, typename SparseLUImpl::GlobalLU_t& glu, SparseLUImpl& luImpl) : m_jcol(jcol), m_jsuper_ref(jsuper), m_glu(glu), m_luImpl(luImpl) {} bool update_segrep(Index /*krep*/, Index /*jj*/) @@ -57,8 +57,8 @@ struct column_dfs_traits : no_assignment_operator Index m_jcol; Index& m_jsuper_ref; - typename SparseLUImpl::GlobalLU_t& m_glu; - SparseLUImpl& m_luImpl; + typename SparseLUImpl::GlobalLU_t& m_glu; + SparseLUImpl& m_luImpl; }; @@ -89,8 +89,10 @@ struct column_dfs_traits : no_assignment_operator * > 0 number of bytes allocated when run out of space * */ -template -Index SparseLUImpl::column_dfs(const Index m, const Index jcol, IndexVector& perm_r, Index maxsuper, Index& nseg, BlockIndexVector lsub_col, IndexVector& segrep, BlockIndexVector repfnz, IndexVector& xprune, IndexVector& marker, IndexVector& parent, IndexVector& xplore, GlobalLU_t& glu) +template +Index SparseLUImpl::column_dfs(const Index m, const Index jcol, IndexVector& perm_r, Index maxsuper, Index& nseg, + BlockIndexVector lsub_col, IndexVector& segrep, BlockIndexVector repfnz, IndexVector& xprune, + IndexVector& marker, IndexVector& parent, IndexVector& xplore, GlobalLU_t& glu) { Index jsuper = glu.supno(jcol); diff --git a/Eigen/src/SparseLU/SparseLU_copy_to_ucol.h b/Eigen/src/SparseLU/SparseLU_copy_to_ucol.h index 170610d9f..bf237951d 100644 --- a/Eigen/src/SparseLU/SparseLU_copy_to_ucol.h +++ b/Eigen/src/SparseLU/SparseLU_copy_to_ucol.h @@ -46,8 +46,9 @@ namespace internal { * > 0 - number of bytes allocated when run out of space * */ -template -Index SparseLUImpl::copy_to_ucol(const Index jcol, const Index nseg, IndexVector& segrep, BlockIndexVector repfnz ,IndexVector& perm_r, BlockScalarVector dense, GlobalLU_t& glu) +template +Index SparseLUImpl::copy_to_ucol(const Index jcol, const Index nseg, IndexVector& segrep, + BlockIndexVector repfnz ,IndexVector& perm_r, BlockScalarVector dense, GlobalLU_t& glu) { Index ksub, krep, ksupno; diff --git a/Eigen/src/SparseLU/SparseLU_gemm_kernel.h b/Eigen/src/SparseLU/SparseLU_gemm_kernel.h index 9e4e3e72b..7420b4d17 100644 --- a/Eigen/src/SparseLU/SparseLU_gemm_kernel.h +++ b/Eigen/src/SparseLU/SparseLU_gemm_kernel.h @@ -21,7 +21,7 @@ namespace internal { * - lda and ldc must be multiples of the respective packet size * - C must have the same alignment as A */ -template +template EIGEN_DONT_INLINE void sparselu_gemm(Index m, Index n, Index d, const Scalar* A, Index lda, const Scalar* B, Index ldb, Scalar* C, Index ldc) { diff --git a/Eigen/src/SparseLU/SparseLU_heap_relax_snode.h b/Eigen/src/SparseLU/SparseLU_heap_relax_snode.h index 7a4e4305a..4092f842f 100644 --- a/Eigen/src/SparseLU/SparseLU_heap_relax_snode.h +++ b/Eigen/src/SparseLU/SparseLU_heap_relax_snode.h @@ -42,8 +42,8 @@ namespace internal { * \param descendants Number of descendants of each node in the etree * \param relax_end last column in a supernode */ -template -void SparseLUImpl::heap_relax_snode (const Index n, IndexVector& et, const Index relax_columns, IndexVector& descendants, IndexVector& relax_end) +template +void SparseLUImpl::heap_relax_snode (const Index n, IndexVector& et, const Index relax_columns, IndexVector& descendants, IndexVector& relax_end) { // The etree may not be postordered, but its heap ordered @@ -75,7 +75,7 @@ void SparseLUImpl::heap_relax_snode (const Index n, IndexVector& e } // Identify the relaxed supernodes by postorder traversal of the etree Index snode_start; // beginning of a snode - Index k; + StorageIndex k; Index nsuper_et_post = 0; // Number of relaxed snodes in postordered etree Index nsuper_et = 0; // Number of relaxed snodes in the original etree Index l; diff --git a/Eigen/src/SparseLU/SparseLU_kernel_bmod.h b/Eigen/src/SparseLU/SparseLU_kernel_bmod.h index cad149ded..9513f8369 100644 --- a/Eigen/src/SparseLU/SparseLU_kernel_bmod.h +++ b/Eigen/src/SparseLU/SparseLU_kernel_bmod.h @@ -30,13 +30,13 @@ namespace internal { */ template struct LU_kernel_bmod { - template + template static EIGEN_DONT_INLINE void run(const Index segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros); }; template -template +template EIGEN_DONT_INLINE void LU_kernel_bmod::run(const Index segsize, BlockScalarVector& dense, ScalarVector& tempv, ScalarVector& lusup, Index& luptr, const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros) { @@ -91,21 +91,22 @@ EIGEN_DONT_INLINE void LU_kernel_bmod::run(const Index seg template <> struct LU_kernel_bmod<1> { - template + template static EIGEN_DONT_INLINE void run(const Index /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr, const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros); }; -template +template EIGEN_DONT_INLINE void LU_kernel_bmod<1>::run(const Index /*segsize*/, BlockScalarVector& dense, ScalarVector& /*tempv*/, ScalarVector& lusup, Index& luptr, const Index lda, const Index nrow, IndexVector& lsub, const Index lptr, const Index no_zeros) { typedef typename ScalarVector::Scalar Scalar; + typedef typename IndexVector::Scalar StorageIndex; Scalar f = dense(lsub(lptr + no_zeros)); luptr += lda * no_zeros + no_zeros + 1; const Scalar* a(lusup.data() + luptr); - const /*typename IndexVector::Scalar*/Index* irow(lsub.data()+lptr + no_zeros + 1); + const StorageIndex* irow(lsub.data()+lptr + no_zeros + 1); Index i = 0; for (; i+1 < nrow; i+=2) { diff --git a/Eigen/src/SparseLU/SparseLU_panel_bmod.h b/Eigen/src/SparseLU/SparseLU_panel_bmod.h index da0e0fc3c..bd3cf87b9 100644 --- a/Eigen/src/SparseLU/SparseLU_panel_bmod.h +++ b/Eigen/src/SparseLU/SparseLU_panel_bmod.h @@ -52,8 +52,8 @@ namespace internal { * * */ -template -void SparseLUImpl::panel_bmod(const Index m, const Index w, const Index jcol, +template +void SparseLUImpl::panel_bmod(const Index m, const Index w, const Index jcol, const Index nseg, ScalarVector& dense, ScalarVector& tempv, IndexVector& segrep, IndexVector& repfnz, GlobalLU_t& glu) { diff --git a/Eigen/src/SparseLU/SparseLU_panel_dfs.h b/Eigen/src/SparseLU/SparseLU_panel_dfs.h index dc0054efd..f4a908ee5 100644 --- a/Eigen/src/SparseLU/SparseLU_panel_dfs.h +++ b/Eigen/src/SparseLU/SparseLU_panel_dfs.h @@ -37,8 +37,8 @@ namespace internal { template struct panel_dfs_traits { - typedef typename IndexVector::Scalar Index; - panel_dfs_traits(Index jcol, Index* marker) + typedef typename IndexVector::Scalar StorageIndex; + panel_dfs_traits(Index jcol, StorageIndex* marker) : m_jcol(jcol), m_marker(marker) {} bool update_segrep(Index krep, Index jj) @@ -53,13 +53,13 @@ struct panel_dfs_traits void mem_expand(IndexVector& /*glu.lsub*/, Index /*nextl*/, Index /*chmark*/) {} enum { ExpandMem = false }; Index m_jcol; - Index* m_marker; + StorageIndex* m_marker; }; -template +template template -void SparseLUImpl::dfs_kernel(const Index jj, IndexVector& perm_r, +void SparseLUImpl::dfs_kernel(const Index jj, IndexVector& perm_r, Index& nseg, IndexVector& panel_lsub, IndexVector& segrep, Ref repfnz_col, IndexVector& xprune, Ref marker, IndexVector& parent, IndexVector& xplore, GlobalLU_t& glu, @@ -215,8 +215,8 @@ void SparseLUImpl::dfs_kernel(const Index jj, IndexVector& perm_r, * */ -template -void SparseLUImpl::panel_dfs(const Index m, const Index w, const Index jcol, MatrixType& A, IndexVector& perm_r, Index& nseg, ScalarVector& dense, IndexVector& panel_lsub, IndexVector& segrep, IndexVector& repfnz, IndexVector& xprune, IndexVector& marker, IndexVector& parent, IndexVector& xplore, GlobalLU_t& glu) +template +void SparseLUImpl::panel_dfs(const Index m, const Index w, const Index jcol, MatrixType& A, IndexVector& perm_r, Index& nseg, ScalarVector& dense, IndexVector& panel_lsub, IndexVector& segrep, IndexVector& repfnz, IndexVector& xprune, IndexVector& marker, IndexVector& parent, IndexVector& xplore, GlobalLU_t& glu) { Index nextl_col; // Next available position in panel_lsub[*,jj] diff --git a/Eigen/src/SparseLU/SparseLU_pivotL.h b/Eigen/src/SparseLU/SparseLU_pivotL.h index 457789c78..01f5ba4e9 100644 --- a/Eigen/src/SparseLU/SparseLU_pivotL.h +++ b/Eigen/src/SparseLU/SparseLU_pivotL.h @@ -56,8 +56,8 @@ namespace internal { * \return 0 if success, i > 0 if U(i,i) is exactly zero * */ -template -Index SparseLUImpl::pivotL(const Index jcol, const RealScalar& diagpivotthresh, IndexVector& perm_r, IndexVector& iperm_c, Index& pivrow, GlobalLU_t& glu) +template +Index SparseLUImpl::pivotL(const Index jcol, const RealScalar& diagpivotthresh, IndexVector& perm_r, IndexVector& iperm_c, Index& pivrow, GlobalLU_t& glu) { Index fsupc = (glu.xsup)((glu.supno)(jcol)); // First column in the supernode containing the column jcol @@ -67,7 +67,7 @@ Index SparseLUImpl::pivotL(const Index jcol, const RealScalar& dia Index lda = glu.xlusup(fsupc+1) - glu.xlusup(fsupc); // leading dimension Scalar* lu_sup_ptr = &(glu.lusup.data()[glu.xlusup(fsupc)]); // Start of the current supernode Scalar* lu_col_ptr = &(glu.lusup.data()[glu.xlusup(jcol)]); // Start of jcol in the supernode - Index* lsub_ptr = &(glu.lsub.data()[lptr]); // Start of row indices of the supernode + StorageIndex* lsub_ptr = &(glu.lsub.data()[lptr]); // Start of row indices of the supernode // Determine the largest abs numerical value for partial pivoting Index diagind = iperm_c(jcol); // diagonal index diff --git a/Eigen/src/SparseLU/SparseLU_pruneL.h b/Eigen/src/SparseLU/SparseLU_pruneL.h index 66460d168..13133fcc2 100644 --- a/Eigen/src/SparseLU/SparseLU_pruneL.h +++ b/Eigen/src/SparseLU/SparseLU_pruneL.h @@ -49,8 +49,9 @@ namespace internal { * \param glu Global LU data * */ -template -void SparseLUImpl::pruneL(const Index jcol, const IndexVector& perm_r, const Index pivrow, const Index nseg, const IndexVector& segrep, BlockIndexVector repfnz, IndexVector& xprune, GlobalLU_t& glu) +template +void SparseLUImpl::pruneL(const Index jcol, const IndexVector& perm_r, const Index pivrow, const Index nseg, + const IndexVector& segrep, BlockIndexVector repfnz, IndexVector& xprune, GlobalLU_t& glu) { // For each supernode-rep irep in U(*,j] Index jsupno = glu.supno(jcol); diff --git a/Eigen/src/SparseLU/SparseLU_relax_snode.h b/Eigen/src/SparseLU/SparseLU_relax_snode.h index 58ec32e27..21c182d56 100644 --- a/Eigen/src/SparseLU/SparseLU_relax_snode.h +++ b/Eigen/src/SparseLU/SparseLU_relax_snode.h @@ -43,8 +43,8 @@ namespace internal { * \param descendants Number of descendants of each node in the etree * \param relax_end last column in a supernode */ -template -void SparseLUImpl::relax_snode (const Index n, IndexVector& et, const Index relax_columns, IndexVector& descendants, IndexVector& relax_end) +template +void SparseLUImpl::relax_snode (const Index n, IndexVector& et, const Index relax_columns, IndexVector& descendants, IndexVector& relax_end) { // compute the number of descendants of each node in the etree -- cgit v1.2.3