aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/OrderingMethods/Ordering.h4
-rw-r--r--Eigen/src/SparseCholesky/SimplicialCholesky.h5
-rw-r--r--test/bicgstab.cpp2
-rw-r--r--test/conjugate_gradient.cpp18
-rw-r--r--test/simplicial_cholesky.cpp24
5 files changed, 29 insertions, 24 deletions
diff --git a/Eigen/src/OrderingMethods/Ordering.h b/Eigen/src/OrderingMethods/Ordering.h
index e88e637a4..cb838d04a 100644
--- a/Eigen/src/OrderingMethods/Ordering.h
+++ b/Eigen/src/OrderingMethods/Ordering.h
@@ -90,11 +90,11 @@ class AMDOrdering
* \note Returns an empty permutation matrix
* \tparam Index The type of indices of the matrix
*/
-template <typename Index>
+template <typename StorageIndex>
class NaturalOrdering
{
public:
- typedef PermutationMatrix<Dynamic, Dynamic, Index> PermutationType;
+ typedef PermutationMatrix<Dynamic, Dynamic, StorageIndex> PermutationType;
/** Compute the permutation vector from a column-major sparse matrix */
template <typename MatrixType>
diff --git a/Eigen/src/SparseCholesky/SimplicialCholesky.h b/Eigen/src/SparseCholesky/SimplicialCholesky.h
index a0815e708..f56298e8c 100644
--- a/Eigen/src/SparseCholesky/SimplicialCholesky.h
+++ b/Eigen/src/SparseCholesky/SimplicialCholesky.h
@@ -69,6 +69,7 @@ class SimplicialCholeskyBase : public SparseSolverBase<Derived>
typedef SparseMatrix<Scalar,ColMajor,StorageIndex> CholMatrixType;
typedef CholMatrixType const * ConstCholMatrixPtr;
typedef Matrix<Scalar,Dynamic,1> VectorType;
+ typedef Matrix<StorageIndex,Dynamic,1> VectorI;
public:
@@ -250,8 +251,8 @@ class SimplicialCholeskyBase : public SparseSolverBase<Derived>
CholMatrixType m_matrix;
VectorType m_diag; // the diagonal coefficients (LDLT mode)
- VectorXi m_parent; // elimination tree
- VectorXi m_nonZerosPerCol;
+ VectorI m_parent; // elimination tree
+ VectorI m_nonZerosPerCol;
PermutationMatrix<Dynamic,Dynamic,StorageIndex> m_P; // the permutation
PermutationMatrix<Dynamic,Dynamic,StorageIndex> m_Pinv; // the inverse permutation
diff --git a/test/bicgstab.cpp b/test/bicgstab.cpp
index 6d76389ce..7a9a11330 100644
--- a/test/bicgstab.cpp
+++ b/test/bicgstab.cpp
@@ -26,6 +26,6 @@ template<typename T, typename I> void test_bicgstab_T()
void test_bicgstab()
{
CALL_SUBTEST_1((test_bicgstab_T<double,int>()) );
- CALL_SUBTEST_1((test_bicgstab_T<double,long int>()));
CALL_SUBTEST_2((test_bicgstab_T<std::complex<double>, int>()));
+ CALL_SUBTEST_3((test_bicgstab_T<double,long int>()));
}
diff --git a/test/conjugate_gradient.cpp b/test/conjugate_gradient.cpp
index 019cc4d64..9622fd86d 100644
--- a/test/conjugate_gradient.cpp
+++ b/test/conjugate_gradient.cpp
@@ -10,13 +10,14 @@
#include "sparse_solver.h"
#include <Eigen/IterativeLinearSolvers>
-template<typename T> void test_conjugate_gradient_T()
+template<typename T, typename I> void test_conjugate_gradient_T()
{
- ConjugateGradient<SparseMatrix<T>, Lower > cg_colmajor_lower_diag;
- ConjugateGradient<SparseMatrix<T>, Upper > cg_colmajor_upper_diag;
- ConjugateGradient<SparseMatrix<T>, Lower|Upper> cg_colmajor_loup_diag;
- ConjugateGradient<SparseMatrix<T>, Lower, IdentityPreconditioner> cg_colmajor_lower_I;
- ConjugateGradient<SparseMatrix<T>, Upper, IdentityPreconditioner> cg_colmajor_upper_I;
+ typedef SparseMatrix<T,0,I> SparseMatrixType;
+ ConjugateGradient<SparseMatrixType, Lower > cg_colmajor_lower_diag;
+ ConjugateGradient<SparseMatrixType, Upper > cg_colmajor_upper_diag;
+ ConjugateGradient<SparseMatrixType, Lower|Upper> cg_colmajor_loup_diag;
+ ConjugateGradient<SparseMatrixType, Lower, IdentityPreconditioner> cg_colmajor_lower_I;
+ ConjugateGradient<SparseMatrixType, Upper, IdentityPreconditioner> cg_colmajor_upper_I;
CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_lower_diag) );
CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_upper_diag) );
@@ -27,6 +28,7 @@ template<typename T> void test_conjugate_gradient_T()
void test_conjugate_gradient()
{
- CALL_SUBTEST_1(test_conjugate_gradient_T<double>());
- CALL_SUBTEST_2(test_conjugate_gradient_T<std::complex<double> >());
+ CALL_SUBTEST_1(( test_conjugate_gradient_T<double,int>() ));
+ CALL_SUBTEST_2(( test_conjugate_gradient_T<std::complex<double>, int>() ));
+ CALL_SUBTEST_3(( test_conjugate_gradient_T<double,long int>() ));
}
diff --git a/test/simplicial_cholesky.cpp b/test/simplicial_cholesky.cpp
index 786468421..b7cc2d351 100644
--- a/test/simplicial_cholesky.cpp
+++ b/test/simplicial_cholesky.cpp
@@ -9,16 +9,17 @@
#include "sparse_solver.h"
-template<typename T> void test_simplicial_cholesky_T()
+template<typename T, typename I> void test_simplicial_cholesky_T()
{
- SimplicialCholesky<SparseMatrix<T>, Lower> chol_colmajor_lower_amd;
- SimplicialCholesky<SparseMatrix<T>, Upper> chol_colmajor_upper_amd;
- SimplicialLLT<SparseMatrix<T>, Lower> llt_colmajor_lower_amd;
- SimplicialLLT<SparseMatrix<T>, Upper> llt_colmajor_upper_amd;
- SimplicialLDLT<SparseMatrix<T>, Lower> ldlt_colmajor_lower_amd;
- SimplicialLDLT<SparseMatrix<T>, Upper> ldlt_colmajor_upper_amd;
- SimplicialLDLT<SparseMatrix<T>, Lower, NaturalOrdering<int> > ldlt_colmajor_lower_nat;
- SimplicialLDLT<SparseMatrix<T>, Upper, NaturalOrdering<int> > ldlt_colmajor_upper_nat;
+ typedef SparseMatrix<T,0,I> SparseMatrixType;
+ SimplicialCholesky<SparseMatrixType, Lower> chol_colmajor_lower_amd;
+ SimplicialCholesky<SparseMatrixType, Upper> chol_colmajor_upper_amd;
+ SimplicialLLT< SparseMatrixType, Lower> llt_colmajor_lower_amd;
+ SimplicialLLT< SparseMatrixType, Upper> llt_colmajor_upper_amd;
+ SimplicialLDLT< SparseMatrixType, Lower> ldlt_colmajor_lower_amd;
+ SimplicialLDLT< SparseMatrixType, Upper> ldlt_colmajor_upper_amd;
+ SimplicialLDLT< SparseMatrixType, Lower, NaturalOrdering<I> > ldlt_colmajor_lower_nat;
+ SimplicialLDLT< SparseMatrixType, Upper, NaturalOrdering<I> > ldlt_colmajor_upper_nat;
check_sparse_spd_solving(chol_colmajor_lower_amd);
check_sparse_spd_solving(chol_colmajor_upper_amd);
@@ -40,6 +41,7 @@ template<typename T> void test_simplicial_cholesky_T()
void test_simplicial_cholesky()
{
- CALL_SUBTEST_1(test_simplicial_cholesky_T<double>());
- CALL_SUBTEST_2(test_simplicial_cholesky_T<std::complex<double> >());
+ CALL_SUBTEST_1(( test_simplicial_cholesky_T<double,int>() ));
+ CALL_SUBTEST_2(( test_simplicial_cholesky_T<std::complex<double>, int>() ));
+ CALL_SUBTEST_3(( test_simplicial_cholesky_T<double,long int>() ));
}