// 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 // // This Source Code Form is subject to the terms of the Mozilla // Public License v. 2.0. If a copy of the MPL was not distributed // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. #include "sparse.h" template void initSPD(double density, Matrix& refMat, SparseMatrix& sparseMat) { Matrix aux(refMat.rows(),refMat.cols()); initSparse(density,refMat,sparseMat); refMat = refMat * refMat.adjoint(); for (int k=0; k<2; ++k) { initSparse(density,aux,sparseMat,ForceNonZeroDiag); refMat += aux * aux.adjoint(); } sparseMat.startFill(); for (int j=0 ; j 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)); // lower - transpose initSparse(density, refMat2, m2, ForceNonZeroDiag|MakeLowerTriangular, &zeroCoords, &nonzeroCoords); VERIFY_IS_APPROX(refMat2.template marked().transpose().solveTriangular(vec2), m2.template marked().transpose().solveTriangular(vec3)); // upper initSparse(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, &zeroCoords, &nonzeroCoords); VERIFY_IS_APPROX(refMat2.template marked().solveTriangular(vec2), m2.template marked().solveTriangular(vec3)); // upper - transpose initSparse(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, &zeroCoords, &nonzeroCoords); VERIFY_IS_APPROX(refMat2.template marked().transpose().solveTriangular(vec2), m2.template marked().transpose().solveTriangular(vec3)); } // test LLT { // 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); initSPD(density, refMat2, m2); refMat2.llt().solve(b, &refX); typedef SparseMatrix SparseSelfAdjointMatrix; if (!NumTraits::IsComplex) { 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 if (!NumTraits::IsComplex) { #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 SparseLDLT::solveInPlace) SparseMatrix m2(rows, cols); DenseMatrix refMat2(rows, cols); DenseVector b = DenseVector::Random(cols); DenseVector refX(cols), x(cols); //initSPD(density, refMat2, m2); initSparse(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, 0, 0); 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); #if defined(EIGEN_SUPERLU_SUPPORT) || defined(EIGEN_UMFPACK_SUPPORT) Scalar refDet = refLu.determinant(); #endif 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_eigen2_sparse_solvers() { for(int i = 0; i < g_repeat; i++) { CALL_SUBTEST_1( sparse_solvers(8, 8) ); CALL_SUBTEST_2( sparse_solvers >(16, 16) ); CALL_SUBTEST_1( sparse_solvers(101, 101) ); } }