svd_fill.h 4.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2014-2015 Gael Guennebaud <gael.guennebaud@inria.fr>
  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. template<typename T>
  10. Array<T,4,1> four_denorms();
  11. template<>
  12. Array4f four_denorms() { return Array4f(5.60844e-39f, -5.60844e-39f, 4.94e-44f, -4.94e-44f); }
  13. template<>
  14. Array4d four_denorms() { return Array4d(5.60844e-313, -5.60844e-313, 4.94e-324, -4.94e-324); }
  15. template<typename T>
  16. Array<T,4,1> four_denorms() { return four_denorms<double>().cast<T>(); }
  17. template<typename MatrixType>
  18. void svd_fill_random(MatrixType &m, int Option = 0)
  19. {
  20. using std::pow;
  21. typedef typename MatrixType::Scalar Scalar;
  22. typedef typename MatrixType::RealScalar RealScalar;
  23. Index diagSize = (std::min)(m.rows(), m.cols());
  24. RealScalar s = std::numeric_limits<RealScalar>::max_exponent10/4;
  25. s = internal::random<RealScalar>(1,s);
  26. Matrix<RealScalar,Dynamic,1> d = Matrix<RealScalar,Dynamic,1>::Random(diagSize);
  27. for(Index k=0; k<diagSize; ++k)
  28. d(k) = d(k)*pow(RealScalar(10),internal::random<RealScalar>(-s,s));
  29. bool dup = internal::random<int>(0,10) < 3;
  30. bool unit_uv = internal::random<int>(0,10) < (dup?7:3); // if we duplicate some diagonal entries, then increase the chance to preserve them using unitary U and V factors
  31. // duplicate some singular values
  32. if(dup)
  33. {
  34. Index n = internal::random<Index>(0,d.size()-1);
  35. for(Index i=0; i<n; ++i)
  36. d(internal::random<Index>(0,d.size()-1)) = d(internal::random<Index>(0,d.size()-1));
  37. }
  38. Matrix<Scalar,Dynamic,Dynamic> U(m.rows(),diagSize);
  39. Matrix<Scalar,Dynamic,Dynamic> VT(diagSize,m.cols());
  40. if(unit_uv)
  41. {
  42. // in very rare cases let's try with a pure diagonal matrix
  43. if(internal::random<int>(0,10) < 1)
  44. {
  45. U.setIdentity();
  46. VT.setIdentity();
  47. }
  48. else
  49. {
  50. createRandomPIMatrixOfRank(diagSize,U.rows(), U.cols(), U);
  51. createRandomPIMatrixOfRank(diagSize,VT.rows(), VT.cols(), VT);
  52. }
  53. }
  54. else
  55. {
  56. U.setRandom();
  57. VT.setRandom();
  58. }
  59. Matrix<Scalar,Dynamic,1> samples(9);
  60. samples << 0, four_denorms<RealScalar>(),
  61. -RealScalar(1)/NumTraits<RealScalar>::highest(), RealScalar(1)/NumTraits<RealScalar>::highest(), (std::numeric_limits<RealScalar>::min)(), pow((std::numeric_limits<RealScalar>::min)(),0.8);
  62. if(Option==Symmetric)
  63. {
  64. m = U * d.asDiagonal() * U.transpose();
  65. // randomly nullify some rows/columns
  66. {
  67. Index count = internal::random<Index>(-diagSize,diagSize);
  68. for(Index k=0; k<count; ++k)
  69. {
  70. Index i = internal::random<Index>(0,diagSize-1);
  71. m.row(i).setZero();
  72. m.col(i).setZero();
  73. }
  74. if(count<0)
  75. // (partly) cancel some coeffs
  76. if(!(dup && unit_uv))
  77. {
  78. Index n = internal::random<Index>(0,m.size()-1);
  79. for(Index k=0; k<n; ++k)
  80. {
  81. Index i = internal::random<Index>(0,m.rows()-1);
  82. Index j = internal::random<Index>(0,m.cols()-1);
  83. m(j,i) = m(i,j) = samples(internal::random<Index>(0,samples.size()-1));
  84. if(NumTraits<Scalar>::IsComplex)
  85. *(&numext::real_ref(m(j,i))+1) = *(&numext::real_ref(m(i,j))+1) = samples.real()(internal::random<Index>(0,samples.size()-1));
  86. }
  87. }
  88. }
  89. }
  90. else
  91. {
  92. m = U * d.asDiagonal() * VT;
  93. // (partly) cancel some coeffs
  94. if(!(dup && unit_uv))
  95. {
  96. Index n = internal::random<Index>(0,m.size()-1);
  97. for(Index k=0; k<n; ++k)
  98. {
  99. Index i = internal::random<Index>(0,m.rows()-1);
  100. Index j = internal::random<Index>(0,m.cols()-1);
  101. m(i,j) = samples(internal::random<Index>(0,samples.size()-1));
  102. if(NumTraits<Scalar>::IsComplex)
  103. *(&numext::real_ref(m(i,j))+1) = samples.real()(internal::random<Index>(0,samples.size()-1));
  104. }
  105. }
  106. }
  107. }