aboutsummaryrefslogtreecommitdiffhomepage
path: root/unsupported/Eigen/src/SparseExtra/Amd.h
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2011-06-06 10:17:28 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2011-06-06 10:17:28 +0200
commit421ece38e1995ec4df12213d6fd567fa18222cca (patch)
tree6966bc7910a6a91f6970e16532a50aab284b8c3e /unsupported/Eigen/src/SparseExtra/Amd.h
parent7a61a564ef7d5403bcf3eef0c84252cc8bf73705 (diff)
Sparse: fix long int as index type in simplicial cholesky and other decompositions
Diffstat (limited to 'unsupported/Eigen/src/SparseExtra/Amd.h')
-rw-r--r--unsupported/Eigen/src/SparseExtra/Amd.h10
1 files changed, 5 insertions, 5 deletions
diff --git a/unsupported/Eigen/src/SparseExtra/Amd.h b/unsupported/Eigen/src/SparseExtra/Amd.h
index 52fd56bc4..3cf8bd1e1 100644
--- a/unsupported/Eigen/src/SparseExtra/Amd.h
+++ b/unsupported/Eigen/src/SparseExtra/Amd.h
@@ -103,7 +103,7 @@ Index cs_tdfs(Index j, Index k, Index *head, const Index *next, Index *post, Ind
* The input matrix \a C must be a selfadjoint compressed column major SparseMatrix object. Both the upper and lower parts have to be stored, but the diagonal entries are optional.
* On exit the values of C are destroyed */
template<typename Scalar, typename Index>
-void minimum_degree_ordering(SparseMatrix<Scalar,ColMajor,Index>& C, PermutationMatrix<Dynamic>& perm)
+void minimum_degree_ordering(SparseMatrix<Scalar,ColMajor,Index>& C, PermutationMatrix<Dynamic,Dynamic,Index>& perm)
{
typedef SparseMatrix<Scalar,ColMajor,Index> CCS;
@@ -151,7 +151,7 @@ void minimum_degree_ordering(SparseMatrix<Scalar,ColMajor,Index>& C, Permutation
elen[i] = 0; // Ek of node i is empty
degree[i] = len[i]; // degree of node i
}
- mark = cs_wclear (0, 0, w, n); /* clear w */
+ mark = cs_wclear<Index>(0, 0, w, n); /* clear w */
elen[n] = -2; /* n is a dead element */
Cp[n] = -1; /* n is a root of assembly tree */
w[n] = 0; /* n is a dead element */
@@ -266,7 +266,7 @@ void minimum_degree_ordering(SparseMatrix<Scalar,ColMajor,Index>& C, Permutation
elen[k] = -2; /* k is now an element */
/* --- Find set differences ----------------------------------------- */
- mark = cs_wclear (mark, lemax, w, n); /* clear w if necessary */
+ mark = cs_wclear<Index>(mark, lemax, w, n); /* clear w if necessary */
for(pk = pk1; pk < pk2; pk++) /* scan 1: find |Le\Lk| */
{
i = Ci[pk];
@@ -349,7 +349,7 @@ void minimum_degree_ordering(SparseMatrix<Scalar,ColMajor,Index>& C, Permutation
} /* scan2 is done */
degree[k] = dk; /* finalize |Lk| */
lemax = std::max<Index>(lemax, dk);
- mark = cs_wclear (mark+lemax, lemax, w, n); /* clear w */
+ mark = cs_wclear<Index>(mark+lemax, lemax, w, n); /* clear w */
/* --- Supernode detection ------------------------------------------ */
for(pk = pk1; pk < pk2; pk++)
@@ -435,7 +435,7 @@ void minimum_degree_ordering(SparseMatrix<Scalar,ColMajor,Index>& C, Permutation
}
for(k = 0, i = 0; i <= n; i++) /* postorder the assembly tree */
{
- if(Cp[i] == -1) k = cs_tdfs (i, k, head, next, perm.indices().data(), w);
+ if(Cp[i] == -1) k = cs_tdfs<Index>(i, k, head, next, perm.indices().data(), w);
}
perm.indices().conservativeResize(n);