jacobisvd.cpp 5.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142
  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. Index rows = m.rows();
  35. Index cols = m.cols();
  36. enum {
  37. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  38. };
  39. MatrixType a = MatrixType::Zero(rows, cols);
  40. a.setZero();
  41. if (ColsAtCompileTime == Dynamic)
  42. {
  43. JacobiSVD<MatrixType, FullPivHouseholderQRPreconditioner> svd_fullqr;
  44. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeFullU|ComputeThinV))
  45. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeThinV))
  46. VERIFY_RAISES_ASSERT(svd_fullqr.compute(a, ComputeThinU|ComputeFullV))
  47. }
  48. }
  49. template<typename MatrixType>
  50. void jacobisvd_method()
  51. {
  52. enum { Size = MatrixType::RowsAtCompileTime };
  53. typedef typename MatrixType::RealScalar RealScalar;
  54. typedef Matrix<RealScalar, Size, 1> RealVecType;
  55. MatrixType m = MatrixType::Identity();
  56. VERIFY_IS_APPROX(m.jacobiSvd().singularValues(), RealVecType::Ones());
  57. VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixU());
  58. VERIFY_RAISES_ASSERT(m.jacobiSvd().matrixV());
  59. VERIFY_IS_APPROX(m.jacobiSvd(ComputeFullU|ComputeFullV).solve(m), m);
  60. }
  61. namespace Foo {
  62. // older compiler require a default constructor for Bar
  63. // cf: https://stackoverflow.com/questions/7411515/
  64. class Bar {public: Bar() {}};
  65. bool operator<(const Bar&, const Bar&) { return true; }
  66. }
  67. // regression test for a very strange MSVC issue for which simply
  68. // including SVDBase.h messes up with std::max and custom scalar type
  69. void msvc_workaround()
  70. {
  71. const Foo::Bar a;
  72. const Foo::Bar b;
  73. std::max EIGEN_NOT_A_MACRO (a,b);
  74. }
  75. void test_jacobisvd()
  76. {
  77. CALL_SUBTEST_3(( jacobisvd_verify_assert(Matrix3f()) ));
  78. CALL_SUBTEST_4(( jacobisvd_verify_assert(Matrix4d()) ));
  79. CALL_SUBTEST_7(( jacobisvd_verify_assert(MatrixXf(10,12)) ));
  80. CALL_SUBTEST_8(( jacobisvd_verify_assert(MatrixXcd(7,5)) ));
  81. CALL_SUBTEST_11(svd_all_trivial_2x2(jacobisvd<Matrix2cd>));
  82. CALL_SUBTEST_12(svd_all_trivial_2x2(jacobisvd<Matrix2d>));
  83. for(int i = 0; i < g_repeat; i++) {
  84. CALL_SUBTEST_3(( jacobisvd<Matrix3f>() ));
  85. CALL_SUBTEST_4(( jacobisvd<Matrix4d>() ));
  86. CALL_SUBTEST_5(( jacobisvd<Matrix<float,3,5> >() ));
  87. CALL_SUBTEST_6(( jacobisvd<Matrix<double,Dynamic,2> >(Matrix<double,Dynamic,2>(10,2)) ));
  88. int r = internal::random<int>(1, 30),
  89. c = internal::random<int>(1, 30);
  90. TEST_SET_BUT_UNUSED_VARIABLE(r)
  91. TEST_SET_BUT_UNUSED_VARIABLE(c)
  92. CALL_SUBTEST_10(( jacobisvd<MatrixXd>(MatrixXd(r,c)) ));
  93. CALL_SUBTEST_7(( jacobisvd<MatrixXf>(MatrixXf(r,c)) ));
  94. CALL_SUBTEST_8(( jacobisvd<MatrixXcd>(MatrixXcd(r,c)) ));
  95. (void) r;
  96. (void) c;
  97. // Test on inf/nan matrix
  98. CALL_SUBTEST_7( (svd_inf_nan<JacobiSVD<MatrixXf>, MatrixXf>()) );
  99. CALL_SUBTEST_10( (svd_inf_nan<JacobiSVD<MatrixXd>, MatrixXd>()) );
  100. // bug1395 test compile-time vectors as input
  101. CALL_SUBTEST_13(( jacobisvd_verify_assert(Matrix<double,6,1>()) ));
  102. CALL_SUBTEST_13(( jacobisvd_verify_assert(Matrix<double,1,6>()) ));
  103. CALL_SUBTEST_13(( jacobisvd_verify_assert(Matrix<double,Dynamic,1>(r)) ));
  104. CALL_SUBTEST_13(( jacobisvd_verify_assert(Matrix<double,1,Dynamic>(c)) ));
  105. }
  106. 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))) ));
  107. 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))) ));
  108. // test matrixbase method
  109. CALL_SUBTEST_1(( jacobisvd_method<Matrix2cd>() ));
  110. CALL_SUBTEST_3(( jacobisvd_method<Matrix3f>() ));
  111. // Test problem size constructors
  112. CALL_SUBTEST_7( JacobiSVD<MatrixXf>(10,10) );
  113. // Check that preallocation avoids subsequent mallocs
  114. CALL_SUBTEST_9( svd_preallocate<void>() );
  115. CALL_SUBTEST_2( svd_underoverflow<void>() );
  116. msvc_workaround();
  117. }