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.

154 lines
4.2 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. #ifndef EIGEN_TESTSPARSE_H
  10. #include "main.h"
  11. #if EIGEN_GNUC_AT_LEAST(4,0) && !defined __ICC
  12. #include <tr1/unordered_map>
  13. #define EIGEN_UNORDERED_MAP_SUPPORT
  14. namespace std {
  15. using std::tr1::unordered_map;
  16. }
  17. #endif
  18. #ifdef EIGEN_GOOGLEHASH_SUPPORT
  19. #include <google/sparse_hash_map>
  20. #endif
  21. #include <Eigen/Cholesky>
  22. #include <Eigen/LU>
  23. #include <Eigen/Sparse>
  24. enum {
  25. ForceNonZeroDiag = 1,
  26. MakeLowerTriangular = 2,
  27. MakeUpperTriangular = 4,
  28. ForceRealDiag = 8
  29. };
  30. /* Initializes both a sparse and dense matrix with same random values,
  31. * and a ratio of \a density non zero entries.
  32. * \param flags is a union of ForceNonZeroDiag, MakeLowerTriangular and MakeUpperTriangular
  33. * allowing to control the shape of the matrix.
  34. * \param zeroCoords and nonzeroCoords allows to get the coordinate lists of the non zero,
  35. * and zero coefficients respectively.
  36. */
  37. template<typename Scalar> void
  38. initSparse(double density,
  39. Matrix<Scalar,Dynamic,Dynamic>& refMat,
  40. SparseMatrix<Scalar>& sparseMat,
  41. int flags = 0,
  42. std::vector<Vector2i>* zeroCoords = 0,
  43. std::vector<Vector2i>* nonzeroCoords = 0)
  44. {
  45. sparseMat.startFill(int(refMat.rows()*refMat.cols()*density));
  46. for(int j=0; j<refMat.cols(); j++)
  47. {
  48. for(int i=0; i<refMat.rows(); i++)
  49. {
  50. Scalar v = (ei_random<double>(0,1) < density) ? ei_random<Scalar>() : Scalar(0);
  51. if ((flags&ForceNonZeroDiag) && (i==j))
  52. {
  53. v = ei_random<Scalar>()*Scalar(3.);
  54. v = v*v + Scalar(5.);
  55. }
  56. if ((flags & MakeLowerTriangular) && j>i)
  57. v = Scalar(0);
  58. else if ((flags & MakeUpperTriangular) && j<i)
  59. v = Scalar(0);
  60. if ((flags&ForceRealDiag) && (i==j))
  61. v = ei_real(v);
  62. if (v!=Scalar(0))
  63. {
  64. sparseMat.fill(i,j) = v;
  65. if (nonzeroCoords)
  66. nonzeroCoords->push_back(Vector2i(i,j));
  67. }
  68. else if (zeroCoords)
  69. {
  70. zeroCoords->push_back(Vector2i(i,j));
  71. }
  72. refMat(i,j) = v;
  73. }
  74. }
  75. sparseMat.endFill();
  76. }
  77. template<typename Scalar> void
  78. initSparse(double density,
  79. Matrix<Scalar,Dynamic,Dynamic>& refMat,
  80. DynamicSparseMatrix<Scalar>& sparseMat,
  81. int flags = 0,
  82. std::vector<Vector2i>* zeroCoords = 0,
  83. std::vector<Vector2i>* nonzeroCoords = 0)
  84. {
  85. sparseMat.startFill(int(refMat.rows()*refMat.cols()*density));
  86. for(int j=0; j<refMat.cols(); j++)
  87. {
  88. for(int i=0; i<refMat.rows(); i++)
  89. {
  90. Scalar v = (ei_random<double>(0,1) < density) ? ei_random<Scalar>() : Scalar(0);
  91. if ((flags&ForceNonZeroDiag) && (i==j))
  92. {
  93. v = ei_random<Scalar>()*Scalar(3.);
  94. v = v*v + Scalar(5.);
  95. }
  96. if ((flags & MakeLowerTriangular) && j>i)
  97. v = Scalar(0);
  98. else if ((flags & MakeUpperTriangular) && j<i)
  99. v = Scalar(0);
  100. if ((flags&ForceRealDiag) && (i==j))
  101. v = ei_real(v);
  102. if (v!=Scalar(0))
  103. {
  104. sparseMat.fill(i,j) = v;
  105. if (nonzeroCoords)
  106. nonzeroCoords->push_back(Vector2i(i,j));
  107. }
  108. else if (zeroCoords)
  109. {
  110. zeroCoords->push_back(Vector2i(i,j));
  111. }
  112. refMat(i,j) = v;
  113. }
  114. }
  115. sparseMat.endFill();
  116. }
  117. template<typename Scalar> void
  118. initSparse(double density,
  119. Matrix<Scalar,Dynamic,1>& refVec,
  120. SparseVector<Scalar>& sparseVec,
  121. std::vector<int>* zeroCoords = 0,
  122. std::vector<int>* nonzeroCoords = 0)
  123. {
  124. sparseVec.reserve(int(refVec.size()*density));
  125. sparseVec.setZero();
  126. for(int i=0; i<refVec.size(); i++)
  127. {
  128. Scalar v = (ei_random<double>(0,1) < density) ? ei_random<Scalar>() : Scalar(0);
  129. if (v!=Scalar(0))
  130. {
  131. sparseVec.fill(i) = v;
  132. if (nonzeroCoords)
  133. nonzeroCoords->push_back(i);
  134. }
  135. else if (zeroCoords)
  136. zeroCoords->push_back(i);
  137. refVec[i] = v;
  138. }
  139. }
  140. #endif // EIGEN_TESTSPARSE_H