You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

200 lines
7.0 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra. Eigen itself is part of the KDE project.
  3. //
  4. // Copyright (C) 2008 Daniel Gomez Ferro <dgomezferro@gmail.com>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "sparse.h"
  10. template<typename Scalar> void
  11. initSPD(double density,
  12. Matrix<Scalar,Dynamic,Dynamic>& refMat,
  13. SparseMatrix<Scalar>& sparseMat)
  14. {
  15. Matrix<Scalar,Dynamic,Dynamic> aux(refMat.rows(),refMat.cols());
  16. initSparse(density,refMat,sparseMat);
  17. refMat = refMat * refMat.adjoint();
  18. for (int k=0; k<2; ++k)
  19. {
  20. initSparse(density,aux,sparseMat,ForceNonZeroDiag);
  21. refMat += aux * aux.adjoint();
  22. }
  23. sparseMat.startFill();
  24. for (int j=0 ; j<sparseMat.cols(); ++j)
  25. for (int i=j ; i<sparseMat.rows(); ++i)
  26. if (refMat(i,j)!=Scalar(0))
  27. sparseMat.fill(i,j) = refMat(i,j);
  28. sparseMat.endFill();
  29. }
  30. template<typename Scalar> void sparse_solvers(int rows, int cols)
  31. {
  32. double density = std::max(8./(rows*cols), 0.01);
  33. typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  34. typedef Matrix<Scalar,Dynamic,1> DenseVector;
  35. // Scalar eps = 1e-6;
  36. DenseVector vec1 = DenseVector::Random(rows);
  37. std::vector<Vector2i> zeroCoords;
  38. std::vector<Vector2i> nonzeroCoords;
  39. // test triangular solver
  40. {
  41. DenseVector vec2 = vec1, vec3 = vec1;
  42. SparseMatrix<Scalar> m2(rows, cols);
  43. DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
  44. // lower
  45. initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag|MakeLowerTriangular, &zeroCoords, &nonzeroCoords);
  46. VERIFY_IS_APPROX(refMat2.template marked<LowerTriangular>().solveTriangular(vec2),
  47. m2.template marked<LowerTriangular>().solveTriangular(vec3));
  48. // lower - transpose
  49. initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag|MakeLowerTriangular, &zeroCoords, &nonzeroCoords);
  50. VERIFY_IS_APPROX(refMat2.template marked<LowerTriangular>().transpose().solveTriangular(vec2),
  51. m2.template marked<LowerTriangular>().transpose().solveTriangular(vec3));
  52. // upper
  53. initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, &zeroCoords, &nonzeroCoords);
  54. VERIFY_IS_APPROX(refMat2.template marked<UpperTriangular>().solveTriangular(vec2),
  55. m2.template marked<UpperTriangular>().solveTriangular(vec3));
  56. // upper - transpose
  57. initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, &zeroCoords, &nonzeroCoords);
  58. VERIFY_IS_APPROX(refMat2.template marked<UpperTriangular>().transpose().solveTriangular(vec2),
  59. m2.template marked<UpperTriangular>().transpose().solveTriangular(vec3));
  60. }
  61. // test LLT
  62. {
  63. // TODO fix the issue with complex (see SparseLLT::solveInPlace)
  64. SparseMatrix<Scalar> m2(rows, cols);
  65. DenseMatrix refMat2(rows, cols);
  66. DenseVector b = DenseVector::Random(cols);
  67. DenseVector refX(cols), x(cols);
  68. initSPD(density, refMat2, m2);
  69. refMat2.llt().solve(b, &refX);
  70. typedef SparseMatrix<Scalar,LowerTriangular|SelfAdjoint> SparseSelfAdjointMatrix;
  71. if (!NumTraits<Scalar>::IsComplex)
  72. {
  73. x = b;
  74. SparseLLT<SparseSelfAdjointMatrix> (m2).solveInPlace(x);
  75. VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LLT: default");
  76. }
  77. #ifdef EIGEN_CHOLMOD_SUPPORT
  78. x = b;
  79. SparseLLT<SparseSelfAdjointMatrix,Cholmod>(m2).solveInPlace(x);
  80. VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LLT: cholmod");
  81. #endif
  82. if (!NumTraits<Scalar>::IsComplex)
  83. {
  84. #ifdef EIGEN_TAUCS_SUPPORT
  85. x = b;
  86. SparseLLT<SparseSelfAdjointMatrix,Taucs>(m2,IncompleteFactorization).solveInPlace(x);
  87. VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LLT: taucs (IncompleteFactorization)");
  88. x = b;
  89. SparseLLT<SparseSelfAdjointMatrix,Taucs>(m2,SupernodalMultifrontal).solveInPlace(x);
  90. VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LLT: taucs (SupernodalMultifrontal)");
  91. x = b;
  92. SparseLLT<SparseSelfAdjointMatrix,Taucs>(m2,SupernodalLeftLooking).solveInPlace(x);
  93. VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LLT: taucs (SupernodalLeftLooking)");
  94. #endif
  95. }
  96. }
  97. // test LDLT
  98. if (!NumTraits<Scalar>::IsComplex)
  99. {
  100. // TODO fix the issue with complex (see SparseLDLT::solveInPlace)
  101. SparseMatrix<Scalar> m2(rows, cols);
  102. DenseMatrix refMat2(rows, cols);
  103. DenseVector b = DenseVector::Random(cols);
  104. DenseVector refX(cols), x(cols);
  105. //initSPD(density, refMat2, m2);
  106. initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, 0, 0);
  107. refMat2 += refMat2.adjoint();
  108. refMat2.diagonal() *= 0.5;
  109. refMat2.ldlt().solve(b, &refX);
  110. typedef SparseMatrix<Scalar,UpperTriangular|SelfAdjoint> SparseSelfAdjointMatrix;
  111. x = b;
  112. SparseLDLT<SparseSelfAdjointMatrix> ldlt(m2);
  113. if (ldlt.succeeded())
  114. ldlt.solveInPlace(x);
  115. VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LDLT: default");
  116. }
  117. // test LU
  118. {
  119. static int count = 0;
  120. SparseMatrix<Scalar> m2(rows, cols);
  121. DenseMatrix refMat2(rows, cols);
  122. DenseVector b = DenseVector::Random(cols);
  123. DenseVector refX(cols), x(cols);
  124. initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag, &zeroCoords, &nonzeroCoords);
  125. LU<DenseMatrix> refLu(refMat2);
  126. refLu.solve(b, &refX);
  127. #if defined(EIGEN_SUPERLU_SUPPORT) || defined(EIGEN_UMFPACK_SUPPORT)
  128. Scalar refDet = refLu.determinant();
  129. #endif
  130. x.setZero();
  131. // // SparseLU<SparseMatrix<Scalar> > (m2).solve(b,&x);
  132. // // VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: default");
  133. #ifdef EIGEN_SUPERLU_SUPPORT
  134. {
  135. x.setZero();
  136. SparseLU<SparseMatrix<Scalar>,SuperLU> slu(m2);
  137. if (slu.succeeded())
  138. {
  139. if (slu.solve(b,&x)) {
  140. VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: SuperLU");
  141. }
  142. // std::cerr << refDet << " == " << slu.determinant() << "\n";
  143. if (count==0) {
  144. VERIFY_IS_APPROX(refDet,slu.determinant()); // FIXME det is not very stable for complex
  145. }
  146. }
  147. }
  148. #endif
  149. #ifdef EIGEN_UMFPACK_SUPPORT
  150. {
  151. // check solve
  152. x.setZero();
  153. SparseLU<SparseMatrix<Scalar>,UmfPack> slu(m2);
  154. if (slu.succeeded()) {
  155. if (slu.solve(b,&x)) {
  156. if (count==0) {
  157. VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: umfpack"); // FIXME solve is not very stable for complex
  158. }
  159. }
  160. VERIFY_IS_APPROX(refDet,slu.determinant());
  161. // TODO check the extracted data
  162. //std::cerr << slu.matrixL() << "\n";
  163. }
  164. }
  165. #endif
  166. count++;
  167. }
  168. }
  169. void test_eigen2_sparse_solvers()
  170. {
  171. for(int i = 0; i < g_repeat; i++) {
  172. CALL_SUBTEST_1( sparse_solvers<double>(8, 8) );
  173. CALL_SUBTEST_2( sparse_solvers<std::complex<double> >(16, 16) );
  174. CALL_SUBTEST_1( sparse_solvers<double>(101, 101) );
  175. }
  176. }