From 86ccd99d8d9a87d03f2f327766a02cc13849b54d Mon Sep 17 00:00:00 2001 From: Gael Guennebaud Date: Wed, 5 Nov 2008 13:47:55 +0000 Subject: Several improvements in sparse module: * add a LDL^T factorization with solver using code from T. Davis's LDL library (LPGL2.1+) * various bug fixes in trianfular solver, matrix product, etc. * improve cmake files for the supported libraries * split the sparse unit test * etc. --- test/sparse_solvers.cpp | 180 ++++++++++++++++++++++++++++++++++++++++++++++++ 1 file changed, 180 insertions(+) create mode 100644 test/sparse_solvers.cpp (limited to 'test/sparse_solvers.cpp') diff --git a/test/sparse_solvers.cpp b/test/sparse_solvers.cpp new file mode 100644 index 000000000..b73494879 --- /dev/null +++ b/test/sparse_solvers.cpp @@ -0,0 +1,180 @@ +// This file is part of Eigen, a lightweight C++ template library +// for linear algebra. Eigen itself is part of the KDE project. +// +// Copyright (C) 2008 Daniel Gomez Ferro +// +// 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 . + +#include "sparse.h" + +template void sparse_solvers(int rows, int cols) +{ + double density = std::max(8./(rows*cols), 0.01); + typedef Matrix DenseMatrix; + typedef Matrix DenseVector; + Scalar eps = 1e-6; + + DenseVector vec1 = DenseVector::Random(rows); + + std::vector zeroCoords; + std::vector nonzeroCoords; + + // test triangular solver + { + DenseVector vec2 = vec1, vec3 = vec1; + SparseMatrix m2(rows, cols); + DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols); + + // lower + initSparse(density, refMat2, m2, ForceNonZeroDiag|MakeLowerTriangular, &zeroCoords, &nonzeroCoords); + VERIFY_IS_APPROX(refMat2.template marked().solveTriangular(vec2), + m2.template marked().solveTriangular(vec3)); + + // upper + initSparse(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, &zeroCoords, &nonzeroCoords); + VERIFY_IS_APPROX(refMat2.template marked().solveTriangular(vec2), + m2.template marked().solveTriangular(vec3)); + + // TODO test row major + } + + // test LLT + if (!NumTraits::IsComplex) + { + // TODO fix the issue with complex (see SparseLLT::solveInPlace) + SparseMatrix m2(rows, cols); + DenseMatrix refMat2(rows, cols); + + DenseVector b = DenseVector::Random(cols); + DenseVector refX(cols), x(cols); + + initSparse(density, refMat2, m2, ForceNonZeroDiag|MakeLowerTriangular, &zeroCoords, &nonzeroCoords); + refMat2 += refMat2.adjoint(); + refMat2.diagonal() *= 0.5; + + refMat2.llt().solve(b, &refX); + typedef SparseMatrix SparseSelfAdjointMatrix; + x = b; + SparseLLT (m2).solveInPlace(x); + //VERIFY(refX.isApprox(x,test_precision()) && "LLT: default"); + #ifdef EIGEN_CHOLMOD_SUPPORT + x = b; + SparseLLT(m2).solveInPlace(x); + VERIFY(refX.isApprox(x,test_precision()) && "LLT: cholmod"); + #endif + #ifdef EIGEN_TAUCS_SUPPORT + x = b; + SparseLLT(m2,IncompleteFactorization).solveInPlace(x); + VERIFY(refX.isApprox(x,test_precision()) && "LLT: taucs (IncompleteFactorization)"); + x = b; + SparseLLT(m2,SupernodalMultifrontal).solveInPlace(x); + VERIFY(refX.isApprox(x,test_precision()) && "LLT: taucs (SupernodalMultifrontal)"); + x = b; + SparseLLT(m2,SupernodalLeftLooking).solveInPlace(x); + VERIFY(refX.isApprox(x,test_precision()) && "LLT: taucs (SupernodalLeftLooking)"); + #endif + } + + // test LDLT + if (!NumTraits::IsComplex) + { + // TODO fix the issue with complex (see SparseLDT::solveInPlace) + SparseMatrix m2(rows, cols); + DenseMatrix refMat2(rows, cols); + + DenseVector b = DenseVector::Random(cols); + DenseVector refX(cols), x(cols); + + initSparse(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, &zeroCoords, &nonzeroCoords); + refMat2 += refMat2.adjoint(); + refMat2.diagonal() *= 0.5; + + refMat2.ldlt().solve(b, &refX); + typedef SparseMatrix SparseSelfAdjointMatrix; + x = b; + SparseLDLT ldlt(m2); + if (ldlt.succeeded()) + ldlt.solveInPlace(x); + VERIFY(refX.isApprox(x,test_precision()) && "LDLT: default"); + } + + // test LU + { + static int count = 0; + SparseMatrix m2(rows, cols); + DenseMatrix refMat2(rows, cols); + + DenseVector b = DenseVector::Random(cols); + DenseVector refX(cols), x(cols); + + initSparse(density, refMat2, m2, ForceNonZeroDiag, &zeroCoords, &nonzeroCoords); + + LU refLu(refMat2); + refLu.solve(b, &refX); + Scalar refDet = refLu.determinant(); + x.setZero(); + // // SparseLU > (m2).solve(b,&x); + // // VERIFY(refX.isApprox(x,test_precision()) && "LU: default"); + #ifdef EIGEN_SUPERLU_SUPPORT + { + x.setZero(); + SparseLU,SuperLU> slu(m2); + if (slu.succeeded()) + { + if (slu.solve(b,&x)) { + VERIFY(refX.isApprox(x,test_precision()) && "LU: SuperLU"); + } + // std::cerr << refDet << " == " << slu.determinant() << "\n"; + if (count==0) { + VERIFY_IS_APPROX(refDet,slu.determinant()); // FIXME det is not very stable for complex + } + } + } + #endif + #ifdef EIGEN_UMFPACK_SUPPORT + { + // check solve + x.setZero(); + SparseLU,UmfPack> slu(m2); + if (slu.succeeded()) { + if (slu.solve(b,&x)) { + if (count==0) { + VERIFY(refX.isApprox(x,test_precision()) && "LU: umfpack"); // FIXME solve is not very stable for complex + } + } + VERIFY_IS_APPROX(refDet,slu.determinant()); + // TODO check the extracted data + //std::cerr << slu.matrixL() << "\n"; + } + } + #endif + count++; + } + +} + +void test_sparse_solvers() +{ + for(int i = 0; i < g_repeat; i++) { + CALL_SUBTEST( sparse_solvers(8, 8) ); + CALL_SUBTEST( sparse_solvers >(16, 16) ); + CALL_SUBTEST( sparse_solvers(33, 33) ); + } +} -- cgit v1.2.3