aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/SparseLU
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2014-10-06 11:42:31 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2014-10-06 11:42:31 +0200
commitfb53ff1edaa96bc1c6090e5ab4982751dbdfc910 (patch)
tree9eacd5b95f0077bc0c6011c76e1b5057e2638e15 /Eigen/src/SparseLU
parent7a1763995396e6119092ce5fc4ca2d536b6acb73 (diff)
Fix SparseLU regarding uncompressed inputs and avoid manual new/delete calls.
Diffstat (limited to 'Eigen/src/SparseLU')
-rw-r--r--Eigen/src/SparseLU/SparseLU.h30
1 files changed, 16 insertions, 14 deletions
diff --git a/Eigen/src/SparseLU/SparseLU.h b/Eigen/src/SparseLU/SparseLU.h
index 1789830b9..b02796293 100644
--- a/Eigen/src/SparseLU/SparseLU.h
+++ b/Eigen/src/SparseLU/SparseLU.h
@@ -364,30 +364,32 @@ void SparseLU<MatrixType, OrderingType>::analyzePattern(const MatrixType& mat)
//TODO It is possible as in SuperLU to compute row and columns scaling vectors to equilibrate the matrix mat.
+ // Firstly, copy the whole input matrix.
+ m_mat = mat;
+
+ // Compute fill-in ordering
OrderingType ord;
- ord(mat,m_perm_c);
+ ord(m_mat,m_perm_c);
// Apply the permutation to the column of the input matrix
- //First copy the whole input matrix.
- m_mat = mat;
- if (m_perm_c.size()) {
+ if (m_perm_c.size())
+ {
m_mat.uncompress(); //NOTE: The effect of this command is only to create the InnerNonzeros pointers. FIXME : This vector is filled but not subsequently used.
- //Then, permute only the column pointers
- const Index * outerIndexPtr;
- if (mat.isCompressed()) outerIndexPtr = mat.outerIndexPtr();
- else
- {
- Index *outerIndexPtr_t = new Index[mat.cols()+1];
- for(Index i = 0; i <= mat.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];
- outerIndexPtr = outerIndexPtr_t;
- }
+ // Then, permute only the column pointers
+ ei_declare_aligned_stack_constructed_variable(Index,outerIndexPtr,mat.cols()+1,mat.isCompressed()?const_cast<Index*>(mat.outerIndexPtr()):0);
+
+ // If the input matrix 'mat' is uncompressed, then the outer-indices do not match the ones of m_mat, and a copy is thus needed.
+ if(!mat.isCompressed())
+ IndexVector::Map(outerIndexPtr, mat.cols()+1) = IndexVector::Map(m_mat.outerIndexPtr(),mat.cols()+1);
+
+ // Apply the permutation and compute the nnz per column.
for (Index i = 0; i < mat.cols(); i++)
{
m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
}
- if(!mat.isCompressed()) delete[] outerIndexPtr;
}
+
// Compute the column elimination tree of the permuted matrix
IndexVector firstRowElt;
internal::coletree(m_mat, m_etree,firstRowElt);