aboutsummaryrefslogtreecommitdiffhomepage
path: root/test
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2011-11-12 14:11:27 +0100
committerGravatar Gael Guennebaud <g.gael@free.fr>2011-11-12 14:11:27 +0100
commit53fa8517245e0136c83b77526b05ce67de232a56 (patch)
tree99dd17062c742eabfc3626a04c38fd6f72e43bc4 /test
parentdcb66d6b403ed2c4341fdb091f2ef22b73ea8b8a (diff)
move sparse solvers from unsupported/ to main Eigen/ and remove the "not stable yet" warning
Diffstat (limited to 'test')
-rw-r--r--test/CMakeLists.txt56
-rw-r--r--test/bicgstab.cpp47
-rw-r--r--test/conjugate_gradient.cpp47
-rw-r--r--test/simplicial_cholesky.cpp57
-rw-r--r--test/sparse_solver.h193
-rw-r--r--test/superlu_support.cpp41
-rw-r--r--test/umfpack_support.cpp41
7 files changed, 481 insertions, 1 deletions
diff --git a/test/CMakeLists.txt b/test/CMakeLists.txt
index 9b4c90f0c..7086f5251 100644
--- a/test/CMakeLists.txt
+++ b/test/CMakeLists.txt
@@ -12,6 +12,46 @@ foreach(i RANGE 1 999)
)
endforeach()
+# configure blas/lapack (use Eigen's ones)
+set(BLAS_FOUND TRUE)
+set(LAPACK_FOUND TRUE)
+set(BLAS_LIBRARIES eigen_blas)
+set(LAPACK_LIBRARIES eigen_lapack)
+
+set(SPARSE_LIBS " ")
+
+find_package(Cholmod)
+if(CHOLMOD_FOUND AND BLAS_FOUND AND LAPACK_FOUND)
+ add_definitions("-DEIGEN_CHOLMOD_SUPPORT")
+ include_directories(${CHOLMOD_INCLUDES})
+ set(SPARSE_LIBS ${SPARSE_LIBS} ${CHOLMOD_LIBRARIES} ${BLAS_LIBRARIES} ${LAPACK_LIBRARIES})
+ ei_add_property(EIGEN_TESTED_BACKENDS "Cholmod, ")
+else()
+ ei_add_property(EIGEN_MISSING_BACKENDS "Cholmod, ")
+endif()
+
+find_package(Umfpack)
+if(UMFPACK_FOUND AND BLAS_FOUND)
+ add_definitions("-DEIGEN_UMFPACK_SUPPORT")
+ include_directories(${UMFPACK_INCLUDES})
+ set(SPARSE_LIBS ${SPARSE_LIBS} ${UMFPACK_LIBRARIES} ${BLAS_LIBRARIES})
+ set(UMFPACK_ALL_LIBS ${UMFPACK_LIBRARIES} ${BLAS_LIBRARIES})
+ ei_add_property(EIGEN_TESTED_BACKENDS "UmfPack, ")
+else()
+ ei_add_property(EIGEN_MISSING_BACKENDS "UmfPack, ")
+endif()
+
+find_package(SuperLU)
+if(SUPERLU_FOUND AND BLAS_FOUND)
+ add_definitions("-DEIGEN_SUPERLU_SUPPORT")
+ include_directories(${SUPERLU_INCLUDES})
+ set(SPARSE_LIBS ${SPARSE_LIBS} ${SUPERLU_LIBRARIES} ${BLAS_LIBRARIES})
+ set(SUPERLU_ALL_LIBS ${SUPERLU_LIBRARIES} ${BLAS_LIBRARIES})
+ ei_add_property(EIGEN_TESTED_BACKENDS "SuperLU, ")
+else()
+ ei_add_property(EIGEN_MISSING_BACKENDS "SuperLU, ")
+endif()
+
find_package(GSL)
if(GSL_FOUND AND GSL_VERSION_MINOR LESS 9)
set(GSL_FOUND "")
@@ -123,7 +163,7 @@ endif(QT4_FOUND)
ei_add_test(sparse_vector)
ei_add_test(sparse_basic)
ei_add_test(sparse_product)
-ei_add_test(sparse_solvers "" "${SPARSE_LIBS}")
+ei_add_test(sparse_solvers)
ei_add_test(umeyama)
ei_add_test(householder)
ei_add_test(swap)
@@ -140,6 +180,20 @@ ei_add_test(sizeoverflow)
ei_add_test(prec_inverse_4x4)
+ei_add_test(simplicial_cholesky)
+ei_add_test(conjugate_gradient)
+ei_add_test(bicgstab)
+
+if(UMFPACK_FOUND)
+ ei_add_test(umfpack_support "" "${UMFPACK_ALL_LIBS}")
+endif()
+
+if(SUPERLU_FOUND)
+ ei_add_test(superlu_support "" "${SUPERLU_ALL_LIBS}")
+endif()
+
+
+
string(TOLOWER "${CMAKE_CXX_COMPILER}" cmake_cxx_compiler_tolower)
if(cmake_cxx_compiler_tolower MATCHES "qcc")
set(CXX_IS_QCC "ON")
diff --git a/test/bicgstab.cpp b/test/bicgstab.cpp
new file mode 100644
index 000000000..222cd6cbd
--- /dev/null
+++ b/test/bicgstab.cpp
@@ -0,0 +1,47 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#include "sparse_solver.h"
+#include <Eigen/IterativeLinearSolvers>
+
+template<typename T> void test_bicgstab_T()
+{
+ BiCGSTAB<SparseMatrix<T>, DiagonalPreconditioner<T> > bicgstab_colmajor_diag;
+ BiCGSTAB<SparseMatrix<T>, IdentityPreconditioner > bicgstab_colmajor_I;
+ //BiCGSTAB<SparseMatrix<T>, IncompleteLU<T> > bicgstab_colmajor_ilu;
+ //BiCGSTAB<SparseMatrix<T>, SSORPreconditioner<T> > bicgstab_colmajor_ssor;
+
+ CALL_SUBTEST( check_sparse_square_solving(bicgstab_colmajor_diag) );
+ CALL_SUBTEST( check_sparse_square_solving(bicgstab_colmajor_I) );
+ //CALL_SUBTEST( check_sparse_square_solving(bicgstab_colmajor_ilu) );
+ //CALL_SUBTEST( check_sparse_square_solving(bicgstab_colmajor_ssor) );
+}
+
+void test_bicgstab()
+{
+ for(int i = 0; i < g_repeat; i++) {
+ CALL_SUBTEST_1(test_bicgstab_T<double>());
+ //CALL_SUBTEST_2(test_bicgstab_T<std::complex<double> >());
+ }
+}
diff --git a/test/conjugate_gradient.cpp b/test/conjugate_gradient.cpp
new file mode 100644
index 000000000..e76327bac
--- /dev/null
+++ b/test/conjugate_gradient.cpp
@@ -0,0 +1,47 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#include "sparse_solver.h"
+#include <Eigen/IterativeLinearSolvers>
+
+template<typename T> 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, IdentityPreconditioner> cg_colmajor_lower_I;
+ ConjugateGradient<SparseMatrix<T>, 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) );
+ CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_lower_I) );
+ CALL_SUBTEST( check_sparse_spd_solving(cg_colmajor_upper_I) );
+}
+
+void test_conjugate_gradient()
+{
+ for(int i = 0; i < g_repeat; i++) {
+ CALL_SUBTEST_1(test_conjugate_gradient_T<double>());
+ CALL_SUBTEST_2(test_conjugate_gradient_T<std::complex<double> >());
+ }
+}
diff --git a/test/simplicial_cholesky.cpp b/test/simplicial_cholesky.cpp
new file mode 100644
index 000000000..b9a55ef36
--- /dev/null
+++ b/test/simplicial_cholesky.cpp
@@ -0,0 +1,57 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#include "sparse_solver.h"
+
+template<typename T> void test_simplicial_cholesky_T()
+{
+ SimplicialCholesky<SparseMatrix<T>, Lower> chol_colmajor_lower;
+ SimplicialCholesky<SparseMatrix<T>, Upper> chol_colmajor_upper;
+ SimplicialLLt<SparseMatrix<T>, Lower> llt_colmajor_lower;
+ SimplicialLDLt<SparseMatrix<T>, Upper> llt_colmajor_upper;
+ SimplicialLDLt<SparseMatrix<T>, Lower> ldlt_colmajor_lower;
+ SimplicialLDLt<SparseMatrix<T>, Upper> ldlt_colmajor_upper;
+
+ check_sparse_spd_solving(chol_colmajor_lower);
+ check_sparse_spd_solving(chol_colmajor_upper);
+ check_sparse_spd_solving(llt_colmajor_lower);
+ check_sparse_spd_solving(llt_colmajor_upper);
+ check_sparse_spd_solving(ldlt_colmajor_lower);
+ check_sparse_spd_solving(ldlt_colmajor_upper);
+
+ check_sparse_spd_determinant(chol_colmajor_lower);
+ check_sparse_spd_determinant(chol_colmajor_upper);
+ check_sparse_spd_determinant(llt_colmajor_lower);
+ check_sparse_spd_determinant(llt_colmajor_upper);
+ check_sparse_spd_determinant(ldlt_colmajor_lower);
+ check_sparse_spd_determinant(ldlt_colmajor_upper);
+}
+
+void test_simplicial_cholesky()
+{
+ for(int i = 0; i < g_repeat; i++) {
+ CALL_SUBTEST_1(test_simplicial_cholesky_T<double>());
+ CALL_SUBTEST_2(test_simplicial_cholesky_T<std::complex<double> >());
+ }
+}
diff --git a/test/sparse_solver.h b/test/sparse_solver.h
new file mode 100644
index 000000000..51bb33a92
--- /dev/null
+++ b/test/sparse_solver.h
@@ -0,0 +1,193 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#include "sparse.h"
+#include <Eigen/SparseCore>
+
+template<typename Solver, typename Rhs, typename DenseMat, typename DenseRhs>
+void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const DenseMat& dA, const DenseRhs& db)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+
+ DenseRhs refX = dA.lu().solve(db);
+
+ Rhs x(b.rows(), b.cols());
+ Rhs oldb = b;
+
+ solver.compute(A);
+ if (solver.info() != Success)
+ {
+ std::cerr << "sparse solver testing: factorization failed (check_sparse_solving)\n";
+ exit(0);
+ return;
+ }
+ x = solver.solve(b);
+ if (solver.info() != Success)
+ {
+ std::cerr << "sparse SPD: solving failed\n";
+ return;
+ }
+ VERIFY(oldb.isApprox(b) && "sparse SPD: the rhs should not be modified!");
+
+ VERIFY(x.isApprox(refX,test_precision<Scalar>()));
+}
+
+template<typename Solver, typename DenseMat>
+void check_sparse_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+ typedef typename Mat::RealScalar RealScalar;
+
+ solver.compute(A);
+ if (solver.info() != Success)
+ {
+ std::cerr << "sparse solver testing: factorization failed (check_sparse_determinant)\n";
+ return;
+ }
+
+ Scalar refDet = dA.determinant();
+ VERIFY_IS_APPROX(refDet,solver.determinant());
+}
+
+
+template<typename Solver, typename DenseMat>
+int generate_sparse_spd_problem(Solver& , typename Solver::MatrixType& A, typename Solver::MatrixType& halfA, DenseMat& dA, int maxSize = 300)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+ typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
+
+ int size = internal::random<int>(1,maxSize);
+ double density = (std::max)(8./(size*size), 0.01);
+
+ Mat M(size, size);
+ DenseMatrix dM(size, size);
+
+ initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
+
+ A = M * M.adjoint();
+ dA = dM * dM.adjoint();
+
+ halfA.resize(size,size);
+ halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
+
+ return size;
+}
+
+template<typename Solver> void check_sparse_spd_solving(Solver& solver)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+ typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
+ typedef Matrix<Scalar,Dynamic,1> DenseVector;
+
+ // generate the problem
+ Mat A, halfA;
+ DenseMatrix dA;
+ int size = generate_sparse_spd_problem(solver, A, halfA, dA);
+
+ // generate the right hand sides
+ int rhsCols = internal::random<int>(1,16);
+ double density = (std::max)(8./(size*rhsCols), 0.1);
+ Mat B(size,rhsCols);
+ DenseVector b = DenseVector::Random(size);
+ DenseMatrix dB(size,rhsCols);
+ initSparse<Scalar>(density, dB, B);
+
+ check_sparse_solving(solver, A, b, dA, b);
+ check_sparse_solving(solver, halfA, b, dA, b);
+ check_sparse_solving(solver, A, dB, dA, dB);
+ check_sparse_solving(solver, halfA, dB, dA, dB);
+ check_sparse_solving(solver, A, B, dA, dB);
+ check_sparse_solving(solver, halfA, B, dA, dB);
+}
+
+template<typename Solver> void check_sparse_spd_determinant(Solver& solver)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+ typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
+
+ // generate the problem
+ Mat A, halfA;
+ DenseMatrix dA;
+ generate_sparse_spd_problem(solver, A, halfA, dA, 30);
+
+ check_sparse_determinant(solver, A, dA);
+ check_sparse_determinant(solver, halfA, dA );
+}
+
+template<typename Solver, typename DenseMat>
+int generate_sparse_square_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+ typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
+
+ int size = internal::random<int>(1,maxSize);
+ double density = (std::max)(8./(size*size), 0.01);
+
+ A.resize(size,size);
+ dA.resize(size,size);
+
+ initSparse<Scalar>(density, dA, A, ForceNonZeroDiag);
+
+ return size;
+}
+
+template<typename Solver> void check_sparse_square_solving(Solver& solver)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+ typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
+ typedef Matrix<Scalar,Dynamic,1> DenseVector;
+
+ int rhsCols = internal::random<int>(1,16);
+
+ Mat A;
+ DenseMatrix dA;
+ int size = generate_sparse_square_problem(solver, A, dA);
+
+ DenseVector b = DenseVector::Random(size);
+ DenseMatrix dB = DenseMatrix::Random(size,rhsCols);
+
+ check_sparse_solving(solver, A, b, dA, b);
+ check_sparse_solving(solver, A, dB, dA, dB);
+}
+
+template<typename Solver> void check_sparse_square_determinant(Solver& solver)
+{
+ typedef typename Solver::MatrixType Mat;
+ typedef typename Mat::Scalar Scalar;
+ typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
+
+ // generate the problem
+ Mat A;
+ DenseMatrix dA;
+ generate_sparse_square_problem(solver, A, dA, 30);
+
+ check_sparse_determinant(solver, A, dA);
+}
diff --git a/test/superlu_support.cpp b/test/superlu_support.cpp
new file mode 100644
index 000000000..e54d72c1c
--- /dev/null
+++ b/test/superlu_support.cpp
@@ -0,0 +1,41 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#include "sparse_solver.h"
+
+#ifdef EIGEN_SUPERLU_SUPPORT
+#include <Eigen/SuperLUSupport>
+#endif
+
+void test_superlu_support()
+{
+ for(int i = 0; i < g_repeat; i++) {
+ SuperLU<SparseMatrix<double> > superlu_double_colmajor;
+ SuperLU<SparseMatrix<std::complex<double> > > superlu_cplxdouble_colmajor;
+ CALL_SUBTEST_1( check_sparse_square_solving(superlu_double_colmajor) );
+ CALL_SUBTEST_2( check_sparse_square_solving(superlu_cplxdouble_colmajor) );
+ CALL_SUBTEST_1( check_sparse_square_determinant(superlu_double_colmajor) );
+ CALL_SUBTEST_2( check_sparse_square_determinant(superlu_cplxdouble_colmajor) );
+ }
+}
diff --git a/test/umfpack_support.cpp b/test/umfpack_support.cpp
new file mode 100644
index 000000000..16f688302
--- /dev/null
+++ b/test/umfpack_support.cpp
@@ -0,0 +1,41 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#include "sparse_solver.h"
+
+#ifdef EIGEN_UMFPACK_SUPPORT
+#include <Eigen/UmfPackSupport>
+#endif
+
+void test_umfpack_support()
+{
+ for(int i = 0; i < g_repeat; i++) {
+ UmfPackLU<SparseMatrix<double> > umfpack_double_colmajor;
+ UmfPackLU<SparseMatrix<std::complex<double> > > umfpack_cplxdouble_colmajor;
+ CALL_SUBTEST_1(check_sparse_square_solving(umfpack_double_colmajor));
+ CALL_SUBTEST_2(check_sparse_square_solving(umfpack_cplxdouble_colmajor));
+ CALL_SUBTEST_1(check_sparse_square_determinant(umfpack_double_colmajor));
+ CALL_SUBTEST_2(check_sparse_square_determinant(umfpack_cplxdouble_colmajor));
+ }
+}