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.

120 lines
4.6 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2014 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. // discard stack allocation as that too bypasses malloc
  11. #define EIGEN_STACK_ALLOCATION_LIMIT 0
  12. #define EIGEN_RUNTIME_NO_MALLOC
  13. #include "main.h"
  14. #include <Eigen/SVD>
  15. #define SVD_DEFAULT(M) JacobiSVD<M>
  16. #define SVD_FOR_MIN_NORM(M) JacobiSVD<M,ColPivHouseholderQRPreconditioner>
  17. #include "svd_common.h"
  18. // Check all variants of JacobiSVD
  19. template<typename MatrixType>
  20. void jacobisvd(const MatrixType& a = MatrixType(), bool pickrandom = true)
  21. {
  22. MatrixType m = a;
  23. if(pickrandom)
  24. svd_fill_random(m);
  25. CALL_SUBTEST(( svd_test_all_computation_options<JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner> >(m, true) )); // check full only
  26. CALL_SUBTEST(( svd_test_all_computation_options<JacobiSVD<MatrixType, ColPivHouseholderQRPreconditioner> >(m, false) ));
  27. CALL_SUBTEST(( svd_test_all_computation_options<JacobiSVD<MatrixType, HouseholderQRPreconditioner> >(m, false) ));
  28. if(m.rows()==m.cols())
  29. CALL_SUBTEST(( svd_test_all_computation_options<JacobiSVD<MatrixType, NoQRPreconditioner> >(m, false) ));
  30. }
  31. template<typename MatrixType> void jacobisvd_verify_assert(const MatrixType& m)
  32. {
  33. svd_verify_assert<JacobiSVD<MatrixType> >(m);
  34. typedef typename MatrixType::Index Index;
  35. Index rows = m.rows();
  36. Index cols = m.cols();
  37. enum {
  38. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  39. };
  40. MatrixType a = MatrixType::Zero(rows, cols);
  41. a.setZero();
  42. if (ColsAtCompileTime == Dynamic)
  43. {
  44. JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner> svd_fullqr;
  45. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeFullU|ComputeThinV))
  46. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeThinV))
  47. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeFullV))
  48. }
  49. }
  50. template<typename MatrixType>
  51. void jacobisvd_method()
  52. {
  53. enum { Size = MatrixType::RowsAtCompileTime };
  54. typedef typename MatrixType::RealScalar RealScalar;
  55. typedef Matrix<RealScalar, Size, 1> RealVecType;
  56. MatrixType m = MatrixType::Identity();
  57. VERIFY_IS_APPROX(m.jacobiSvd().singularValues(), RealVecType::Ones());
  58. VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixU());
  59. VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixV());
  60. VERIFY_IS_APPROX(m.jacobiSvd(ComputeFullU|ComputeFullV).solve(m), m);
  61. }
  62. void test_jacobisvd()
  63. {
  64. CALL_SUBTEST_3(( jacobisvd_verify_assert(Matrix3f()) ));
  65. CALL_SUBTEST_4(( jacobisvd_verify_assert(Matrix4d()) ));
  66. CALL_SUBTEST_7(( jacobisvd_verify_assert(MatrixXf(10,12)) ));
  67. CALL_SUBTEST_8(( jacobisvd_verify_assert(MatrixXcd(7,5)) ));
  68. CALL_SUBTEST_11(svd_all_trivial_2x2(jacobisvd<Matrix2cd>));
  69. CALL_SUBTEST_12(svd_all_trivial_2x2(jacobisvd<Matrix2d>));
  70. for(int i = 0; i < g_repeat; i++) {
  71. CALL_SUBTEST_3(( jacobisvd<Matrix3f>() ));
  72. CALL_SUBTEST_4(( jacobisvd<Matrix4d>() ));
  73. CALL_SUBTEST_5(( jacobisvd<Matrix<float,3,5> >() ));
  74. CALL_SUBTEST_6(( jacobisvd<Matrix<double,Dynamic,2> >(Matrix<double,Dynamic,2>(10,2)) ));
  75. int r = internal::random<int>(1, 30),
  76. c = internal::random<int>(1, 30);
  77. TEST_SET_BUT_UNUSED_VARIABLE(r)
  78. TEST_SET_BUT_UNUSED_VARIABLE(c)
  79. CALL_SUBTEST_10(( jacobisvd<MatrixXd>(MatrixXd(r,c)) ));
  80. CALL_SUBTEST_7(( jacobisvd<MatrixXf>(MatrixXf(r,c)) ));
  81. CALL_SUBTEST_8(( jacobisvd<MatrixXcd>(MatrixXcd(r,c)) ));
  82. (void) r;
  83. (void) c;
  84. // Test on inf/nan matrix
  85. CALL_SUBTEST_7( (svd_inf_nan<JacobiSVD<MatrixXf>, MatrixXf>()) );
  86. CALL_SUBTEST_10( (svd_inf_nan<JacobiSVD<MatrixXd>, MatrixXd>()) );
  87. }
  88. 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))) ));
  89. 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))) ));
  90. // test matrixbase method
  91. CALL_SUBTEST_1(( jacobisvd_method<Matrix2cd>() ));
  92. CALL_SUBTEST_3(( jacobisvd_method<Matrix3f>() ));
  93. // Test problem size constructors
  94. CALL_SUBTEST_7( JacobiSVD<MatrixXf>(10,10) );
  95. // Check that preallocation avoids subsequent mallocs
  96. CALL_SUBTEST_9( svd_preallocate<void>() );
  97. CALL_SUBTEST_2( svd_underoverflow<void>() );
  98. }