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.

198 lines
6.3 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #include "svd_common.h"
  11. template<typename MatrixType, int QRPreconditioner>
  12. void jacobisvd_check_full(const MatrixType& m, const JacobiSVD<MatrixType, QRPreconditioner>& svd)
  13. {
  14. svd_check_full<MatrixType, JacobiSVD<MatrixType, QRPreconditioner > >(m, svd);
  15. }
  16. template<typename MatrixType, int QRPreconditioner>
  17. void jacobisvd_compare_to_full(const MatrixType& m,
  18. unsigned int computationOptions,
  19. const JacobiSVD<MatrixType, QRPreconditioner>& referenceSvd)
  20. {
  21. svd_compare_to_full<MatrixType, JacobiSVD<MatrixType, QRPreconditioner> >(m, computationOptions, referenceSvd);
  22. }
  23. template<typename MatrixType, int QRPreconditioner>
  24. void jacobisvd_solve(const MatrixType& m, unsigned int computationOptions)
  25. {
  26. svd_solve< MatrixType, JacobiSVD< MatrixType, QRPreconditioner > >(m, computationOptions);
  27. }
  28. template<typename MatrixType, int QRPreconditioner>
  29. void jacobisvd_test_all_computation_options(const MatrixType& m)
  30. {
  31. if (QRPreconditioner == NoQRPreconditioner && m.rows() != m.cols())
  32. return;
  33. JacobiSVD< MatrixType, QRPreconditioner > fullSvd(m, ComputeFullU|ComputeFullV);
  34. svd_test_computation_options_1< MatrixType, JacobiSVD< MatrixType, QRPreconditioner > >(m, fullSvd);
  35. if(QRPreconditioner == FullPivHouseholderQRPreconditioner)
  36. return;
  37. svd_test_computation_options_2< MatrixType, JacobiSVD< MatrixType, QRPreconditioner > >(m, fullSvd);
  38. }
  39. template<typename MatrixType>
  40. void jacobisvd(const MatrixType& a = MatrixType(), bool pickrandom = true)
  41. {
  42. MatrixType m = pickrandom ? MatrixType::Random(a.rows(), a.cols()) : a;
  43. jacobisvd_test_all_computation_options<MatrixType, FullPivHouseholderQRPreconditioner>(m);
  44. jacobisvd_test_all_computation_options<MatrixType, ColPivHouseholderQRPreconditioner>(m);
  45. jacobisvd_test_all_computation_options<MatrixType, HouseholderQRPreconditioner>(m);
  46. jacobisvd_test_all_computation_options<MatrixType, NoQRPreconditioner>(m);
  47. }
  48. template<typename MatrixType>
  49. void jacobisvd_verify_assert(const MatrixType& m)
  50. {
  51. svd_verify_assert<MatrixType, JacobiSVD< MatrixType > >(m);
  52. typedef typename MatrixType::Index Index;
  53. Index rows = m.rows();
  54. Index cols = m.cols();
  55. enum {
  56. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  57. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  58. };
  59. MatrixType a = MatrixType::Zero(rows, cols);
  60. a.setZero();
  61. if (ColsAtCompileTime == Dynamic)
  62. {
  63. JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner> svd_fullqr;
  64. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeFullU|ComputeThinV))
  65. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeThinV))
  66. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeFullV))
  67. }
  68. }
  69. template<typename MatrixType>
  70. void jacobisvd_method()
  71. {
  72. enum { Size = MatrixType::RowsAtCompileTime };
  73. typedef typename MatrixType::RealScalar RealScalar;
  74. typedef Matrix<RealScalar, Size, 1> RealVecType;
  75. MatrixType m = MatrixType::Identity();
  76. VERIFY_IS_APPROX(m.jacobiSvd().singularValues(), RealVecType::Ones());
  77. VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixU());
  78. VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixV());
  79. VERIFY_IS_APPROX(m.jacobiSvd(ComputeFullU|ComputeFullV).solve(m), m);
  80. }
  81. template<typename MatrixType>
  82. void jacobisvd_inf_nan()
  83. {
  84. svd_inf_nan<MatrixType, JacobiSVD< MatrixType > >();
  85. }
  86. // Regression test for bug 286: JacobiSVD loops indefinitely with some
  87. // matrices containing denormal numbers.
  88. void jacobisvd_bug286()
  89. {
  90. #if defined __INTEL_COMPILER
  91. // shut up warning #239: floating point underflow
  92. #pragma warning push
  93. #pragma warning disable 239
  94. #endif
  95. Matrix2d M;
  96. M << -7.90884e-313, -4.94e-324,
  97. 0, 5.60844e-313;
  98. #if defined __INTEL_COMPILER
  99. #pragma warning pop
  100. #endif
  101. JacobiSVD<Matrix2d> svd;
  102. svd.compute(M); // just check we don't loop indefinitely
  103. }
  104. void jacobisvd_preallocate()
  105. {
  106. svd_preallocate< JacobiSVD <MatrixXf> >();
  107. }
  108. void test_jacobisvd()
  109. {
  110. CALL_SUBTEST_11(( jacobisvd<Matrix<double,Dynamic,Dynamic> >
  111. (Matrix<double,Dynamic,Dynamic>(16, 6)) ));
  112. CALL_SUBTEST_3(( jacobisvd_verify_assert(Matrix3f()) ));
  113. CALL_SUBTEST_4(( jacobisvd_verify_assert(Matrix4d()) ));
  114. CALL_SUBTEST_7(( jacobisvd_verify_assert(MatrixXf(10,12)) ));
  115. CALL_SUBTEST_8(( jacobisvd_verify_assert(MatrixXcd(7,5)) ));
  116. for(int i = 0; i < g_repeat; i++) {
  117. Matrix2cd m;
  118. m << 0, 1,
  119. 0, 1;
  120. CALL_SUBTEST_1(( jacobisvd(m, false) ));
  121. m << 1, 0,
  122. 1, 0;
  123. CALL_SUBTEST_1(( jacobisvd(m, false) ));
  124. Matrix2d n;
  125. n << 0, 0,
  126. 0, 0;
  127. CALL_SUBTEST_2(( jacobisvd(n, false) ));
  128. n << 0, 0,
  129. 0, 1;
  130. CALL_SUBTEST_2(( jacobisvd(n, false) ));
  131. CALL_SUBTEST_3(( jacobisvd<Matrix3f>() ));
  132. CALL_SUBTEST_4(( jacobisvd<Matrix4d>() ));
  133. CALL_SUBTEST_5(( jacobisvd<Matrix<float,3,5> >() ));
  134. CALL_SUBTEST_6(( jacobisvd<Matrix<double,Dynamic,2> >(Matrix<double,Dynamic,2>(10,2)) ));
  135. int r = internal::random<int>(1, 30),
  136. c = internal::random<int>(1, 30);
  137. CALL_SUBTEST_7(( jacobisvd<MatrixXf>(MatrixXf(r,c)) ));
  138. CALL_SUBTEST_8(( jacobisvd<MatrixXcd>(MatrixXcd(r,c)) ));
  139. (void) r;
  140. (void) c;
  141. // Test on inf/nan matrix
  142. CALL_SUBTEST_7( jacobisvd_inf_nan<MatrixXf>() );
  143. }
  144. CALL_SUBTEST_7(( jacobisvd<MatrixXf>(MatrixXf(internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/2), internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/2))) ));
  145. CALL_SUBTEST_8(( jacobisvd<MatrixXcd>(MatrixXcd(internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/3), internal::random<int>(EIGEN_TEST_MAX_SIZE/4, EIGEN_TEST_MAX_SIZE/3))) ));
  146. // test matrixbase method
  147. CALL_SUBTEST_1(( jacobisvd_method<Matrix2cd>() ));
  148. CALL_SUBTEST_3(( jacobisvd_method<Matrix3f>() ));
  149. // Test problem size constructors
  150. CALL_SUBTEST_7( JacobiSVD<MatrixXf>(10,10) );
  151. // Check that preallocation avoids subsequent mallocs
  152. CALL_SUBTEST_9( jacobisvd_preallocate() );
  153. // Regression check for bug 286
  154. CALL_SUBTEST_2( jacobisvd_bug286() );
  155. }