bdcsvd.cpp 4.3 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2013 Gauthier Brun <brun.gauthier@gmail.com>
  5. // Copyright (C) 2013 Nicolas Carre <nicolas.carre@ensimag.fr>
  6. // Copyright (C) 2013 Jean Ceccato <jean.ceccato@ensimag.fr>
  7. // Copyright (C) 2013 Pierre Zoppitelli <pierre.zoppitelli@ensimag.fr>
  8. //
  9. // This Source Code Form is subject to the terms of the Mozilla
  10. // Public License v. 2.0. If a copy of the MPL was not distributed
  11. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/
  12. // discard stack allocation as that too bypasses malloc
  13. #define EIGEN_STACK_ALLOCATION_LIMIT 0
  14. #define EIGEN_RUNTIME_NO_MALLOC
  15. #include "main.h"
  16. #include <Eigen/SVD>
  17. #include <iostream>
  18. #include <Eigen/LU>
  19. #define SVD_DEFAULT(M) BDCSVD<M>
  20. #define SVD_FOR_MIN_NORM(M) BDCSVD<M>
  21. #include "svd_common.h"
  22. // Check all variants of JacobiSVD
  23. template<typename MatrixType>
  24. void bdcsvd(const MatrixType& a = MatrixType(), bool pickrandom = true)
  25. {
  26. MatrixType m;
  27. if(pickrandom) {
  28. m.resizeLike(a);
  29. svd_fill_random(m);
  30. }
  31. else
  32. m = a;
  33. CALL_SUBTEST(( svd_test_all_computation_options<BDCSVD<MatrixType> >(m, false) ));
  34. }
  35. template<typename MatrixType>
  36. void bdcsvd_method()
  37. {
  38. enum { Size = MatrixType::RowsAtCompileTime };
  39. typedef typename MatrixType::RealScalar RealScalar;
  40. typedef Matrix<RealScalar, Size, 1> RealVecType;
  41. MatrixType m = MatrixType::Identity();
  42. VERIFY_IS_APPROX(m.bdcSvd().singularValues(), RealVecType::Ones());
  43. VERIFY_RAISES_ASSERT(m.bdcSvd().matrixU());
  44. VERIFY_RAISES_ASSERT(m.bdcSvd().matrixV());
  45. VERIFY_IS_APPROX(m.bdcSvd(ComputeFullU|ComputeFullV).solve(m), m);
  46. VERIFY_IS_APPROX(m.bdcSvd(ComputeFullU|ComputeFullV).transpose().solve(m), m);
  47. VERIFY_IS_APPROX(m.bdcSvd(ComputeFullU|ComputeFullV).adjoint().solve(m), m);
  48. }
  49. // compare the Singular values returned with Jacobi and Bdc
  50. template<typename MatrixType>
  51. void compare_bdc_jacobi(const MatrixType& a = MatrixType(), unsigned int computationOptions = 0)
  52. {
  53. MatrixType m = MatrixType::Random(a.rows(), a.cols());
  54. BDCSVD<MatrixType> bdc_svd(m);
  55. JacobiSVD<MatrixType> jacobi_svd(m);
  56. VERIFY_IS_APPROX(bdc_svd.singularValues(), jacobi_svd.singularValues());
  57. if(computationOptions & ComputeFullU) VERIFY_IS_APPROX(bdc_svd.matrixU(), jacobi_svd.matrixU());
  58. if(computationOptions & ComputeThinU) VERIFY_IS_APPROX(bdc_svd.matrixU(), jacobi_svd.matrixU());
  59. if(computationOptions & ComputeFullV) VERIFY_IS_APPROX(bdc_svd.matrixV(), jacobi_svd.matrixV());
  60. if(computationOptions & ComputeThinV) VERIFY_IS_APPROX(bdc_svd.matrixV(), jacobi_svd.matrixV());
  61. }
  62. EIGEN_DECLARE_TEST(bdcsvd)
  63. {
  64. CALL_SUBTEST_3(( svd_verify_assert<BDCSVD<Matrix3f> >(Matrix3f()) ));
  65. CALL_SUBTEST_4(( svd_verify_assert<BDCSVD<Matrix4d> >(Matrix4d()) ));
  66. CALL_SUBTEST_7(( svd_verify_assert<BDCSVD<MatrixXf> >(MatrixXf(10,12)) ));
  67. CALL_SUBTEST_8(( svd_verify_assert<BDCSVD<MatrixXcd> >(MatrixXcd(7,5)) ));
  68. CALL_SUBTEST_101(( svd_all_trivial_2x2(bdcsvd<Matrix2cd>) ));
  69. CALL_SUBTEST_102(( svd_all_trivial_2x2(bdcsvd<Matrix2d>) ));
  70. for(int i = 0; i < g_repeat; i++) {
  71. CALL_SUBTEST_3(( bdcsvd<Matrix3f>() ));
  72. CALL_SUBTEST_4(( bdcsvd<Matrix4d>() ));
  73. CALL_SUBTEST_5(( bdcsvd<Matrix<float,3,5> >() ));
  74. int r = internal::random<int>(1, EIGEN_TEST_MAX_SIZE/2),
  75. c = internal::random<int>(1, EIGEN_TEST_MAX_SIZE/2);
  76. TEST_SET_BUT_UNUSED_VARIABLE(r)
  77. TEST_SET_BUT_UNUSED_VARIABLE(c)
  78. CALL_SUBTEST_6(( bdcsvd(Matrix<double,Dynamic,2>(r,2)) ));
  79. CALL_SUBTEST_7(( bdcsvd(MatrixXf(r,c)) ));
  80. CALL_SUBTEST_7(( compare_bdc_jacobi(MatrixXf(r,c)) ));
  81. CALL_SUBTEST_10(( bdcsvd(MatrixXd(r,c)) ));
  82. CALL_SUBTEST_10(( compare_bdc_jacobi(MatrixXd(r,c)) ));
  83. CALL_SUBTEST_8(( bdcsvd(MatrixXcd(r,c)) ));
  84. CALL_SUBTEST_8(( compare_bdc_jacobi(MatrixXcd(r,c)) ));
  85. // Test on inf/nan matrix
  86. CALL_SUBTEST_7( (svd_inf_nan<BDCSVD<MatrixXf>, MatrixXf>()) );
  87. CALL_SUBTEST_10( (svd_inf_nan<BDCSVD<MatrixXd>, MatrixXd>()) );
  88. }
  89. // test matrixbase method
  90. CALL_SUBTEST_1(( bdcsvd_method<Matrix2cd>() ));
  91. CALL_SUBTEST_3(( bdcsvd_method<Matrix3f>() ));
  92. // Test problem size constructors
  93. CALL_SUBTEST_7( BDCSVD<MatrixXf>(10,10) );
  94. // Check that preallocation avoids subsequent mallocs
  95. // Disabled because not supported by BDCSVD
  96. // CALL_SUBTEST_9( svd_preallocate<void>() );
  97. CALL_SUBTEST_2( svd_underoverflow<void>() );
  98. }