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.

261 lines
8.7 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. // Copyright (C) 2013 Gauthier Brun <brun.gauthier@gmail.com>
  8. // Copyright (C) 2013 Nicolas Carre <nicolas.carre@ensimag.fr>
  9. // Copyright (C) 2013 Jean Ceccato <jean.ceccato@ensimag.fr>
  10. // Copyright (C) 2013 Pierre Zoppitelli <pierre.zoppitelli@ensimag.fr>
  11. //
  12. // This Source Code Form is subject to the terms of the Mozilla
  13. // Public License v. 2.0. If a copy of the MPL was not distributed
  14. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  15. // discard stack allocation as that too bypasses malloc
  16. #define EIGEN_STACK_ALLOCATION_LIMIT 0
  17. #define EIGEN_RUNTIME_NO_MALLOC
  18. #include "main.h"
  19. #include <unsupported/Eigen/SVD>
  20. #include <Eigen/LU>
  21. // check if "svd" is the good image of "m"
  22. template<typename MatrixType, typename SVD>
  23. void svd_check_full(const MatrixType& m, const SVD& svd)
  24. {
  25. typedef typename MatrixType::Index Index;
  26. Index rows = m.rows();
  27. Index cols = m.cols();
  28. enum {
  29. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  30. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  31. };
  32. typedef typename MatrixType::Scalar Scalar;
  33. typedef Matrix<Scalar, RowsAtCompileTime, RowsAtCompileTime> MatrixUType;
  34. typedef Matrix<Scalar, ColsAtCompileTime, ColsAtCompileTime> MatrixVType;
  35. MatrixType sigma = MatrixType::Zero(rows, cols);
  36. sigma.diagonal() = svd.singularValues().template cast<Scalar>();
  37. MatrixUType u = svd.matrixU();
  38. MatrixVType v = svd.matrixV();
  39. VERIFY_IS_APPROX(m, u * sigma * v.adjoint());
  40. VERIFY_IS_UNITARY(u);
  41. VERIFY_IS_UNITARY(v);
  42. } // end svd_check_full
  43. // Compare to a reference value
  44. template<typename MatrixType, typename SVD>
  45. void svd_compare_to_full(const MatrixType& m,
  46. unsigned int computationOptions,
  47. const SVD& referenceSvd)
  48. {
  49. typedef typename MatrixType::Index Index;
  50. Index rows = m.rows();
  51. Index cols = m.cols();
  52. Index diagSize = (std::min)(rows, cols);
  53. SVD svd(m, computationOptions);
  54. VERIFY_IS_APPROX(svd.singularValues(), referenceSvd.singularValues());
  55. if(computationOptions & ComputeFullU)
  56. VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU());
  57. if(computationOptions & ComputeThinU)
  58. VERIFY_IS_APPROX(svd.matrixU(), referenceSvd.matrixU().leftCols(diagSize));
  59. if(computationOptions & ComputeFullV)
  60. VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV());
  61. if(computationOptions & ComputeThinV)
  62. VERIFY_IS_APPROX(svd.matrixV(), referenceSvd.matrixV().leftCols(diagSize));
  63. } // end svd_compare_to_full
  64. template<typename MatrixType, typename SVD>
  65. void svd_solve(const MatrixType& m, unsigned int computationOptions)
  66. {
  67. typedef typename MatrixType::Scalar Scalar;
  68. typedef typename MatrixType::Index Index;
  69. Index rows = m.rows();
  70. Index cols = m.cols();
  71. enum {
  72. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  73. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  74. };
  75. typedef Matrix<Scalar, RowsAtCompileTime, Dynamic> RhsType;
  76. typedef Matrix<Scalar, ColsAtCompileTime, Dynamic> SolutionType;
  77. RhsType rhs = RhsType::Random(rows, internal::random<Index>(1, cols));
  78. SVD svd(m, computationOptions);
  79. SolutionType x = svd.solve(rhs);
  80. // evaluate normal equation which works also for least-squares solutions
  81. VERIFY_IS_APPROX(m.adjoint()*m*x,m.adjoint()*rhs);
  82. } // end svd_solve
  83. // test computations options
  84. // 2 functions because Jacobisvd can return before the second function
  85. template<typename MatrixType, typename SVD>
  86. void svd_test_computation_options_1(const MatrixType& m, const SVD& fullSvd)
  87. {
  88. svd_check_full< MatrixType, SVD >(m, fullSvd);
  89. svd_solve< MatrixType, SVD >(m, ComputeFullU | ComputeFullV);
  90. }
  91. template<typename MatrixType, typename SVD>
  92. void svd_test_computation_options_2(const MatrixType& m, const SVD& fullSvd)
  93. {
  94. svd_compare_to_full< MatrixType, SVD >(m, ComputeFullU, fullSvd);
  95. svd_compare_to_full< MatrixType, SVD >(m, ComputeFullV, fullSvd);
  96. svd_compare_to_full< MatrixType, SVD >(m, 0, fullSvd);
  97. if (MatrixType::ColsAtCompileTime == Dynamic) {
  98. // thin U/V are only available with dynamic number of columns
  99. svd_compare_to_full< MatrixType, SVD >(m, ComputeFullU|ComputeThinV, fullSvd);
  100. svd_compare_to_full< MatrixType, SVD >(m, ComputeThinV, fullSvd);
  101. svd_compare_to_full< MatrixType, SVD >(m, ComputeThinU|ComputeFullV, fullSvd);
  102. svd_compare_to_full< MatrixType, SVD >(m, ComputeThinU , fullSvd);
  103. svd_compare_to_full< MatrixType, SVD >(m, ComputeThinU|ComputeThinV, fullSvd);
  104. svd_solve<MatrixType, SVD>(m, ComputeFullU | ComputeThinV);
  105. svd_solve<MatrixType, SVD>(m, ComputeThinU | ComputeFullV);
  106. svd_solve<MatrixType, SVD>(m, ComputeThinU | ComputeThinV);
  107. typedef typename MatrixType::Index Index;
  108. Index diagSize = (std::min)(m.rows(), m.cols());
  109. SVD svd(m, ComputeThinU | ComputeThinV);
  110. VERIFY_IS_APPROX(m, svd.matrixU().leftCols(diagSize) * svd.singularValues().asDiagonal() * svd.matrixV().leftCols(diagSize).adjoint());
  111. }
  112. }
  113. template<typename MatrixType, typename SVD>
  114. void svd_verify_assert(const MatrixType& m)
  115. {
  116. typedef typename MatrixType::Scalar Scalar;
  117. typedef typename MatrixType::Index Index;
  118. Index rows = m.rows();
  119. Index cols = m.cols();
  120. enum {
  121. RowsAtCompileTime = MatrixType::RowsAtCompileTime,
  122. ColsAtCompileTime = MatrixType::ColsAtCompileTime
  123. };
  124. typedef Matrix<Scalar, RowsAtCompileTime, 1> RhsType;
  125. RhsType rhs(rows);
  126. SVD svd;
  127. VERIFY_RAISES_ASSERT(svd.matrixU())
  128. VERIFY_RAISES_ASSERT(svd.singularValues())
  129. VERIFY_RAISES_ASSERT(svd.matrixV())
  130. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  131. MatrixType a = MatrixType::Zero(rows, cols);
  132. a.setZero();
  133. svd.compute(a, 0);
  134. VERIFY_RAISES_ASSERT(svd.matrixU())
  135. VERIFY_RAISES_ASSERT(svd.matrixV())
  136. svd.singularValues();
  137. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  138. if (ColsAtCompileTime == Dynamic)
  139. {
  140. svd.compute(a, ComputeThinU);
  141. svd.matrixU();
  142. VERIFY_RAISES_ASSERT(svd.matrixV())
  143. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  144. svd.compute(a, ComputeThinV);
  145. svd.matrixV();
  146. VERIFY_RAISES_ASSERT(svd.matrixU())
  147. VERIFY_RAISES_ASSERT(svd.solve(rhs))
  148. }
  149. else
  150. {
  151. VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinU))
  152. VERIFY_RAISES_ASSERT(svd.compute(a, ComputeThinV))
  153. }
  154. }
  155. // work around stupid msvc error when constructing at compile time an expression that involves
  156. // a division by zero, even if the numeric type has floating point
  157. template<typename Scalar>
  158. EIGEN_DONT_INLINE Scalar zero() { return Scalar(0); }
  159. // workaround aggressive optimization in ICC
  160. template<typename T> EIGEN_DONT_INLINE T sub(T a, T b) { return a - b; }
  161. template<typename MatrixType, typename SVD>
  162. void svd_inf_nan()
  163. {
  164. // all this function does is verify we don't iterate infinitely on nan/inf values
  165. SVD svd;
  166. typedef typename MatrixType::Scalar Scalar;
  167. Scalar some_inf = Scalar(1) / zero<Scalar>();
  168. VERIFY(sub(some_inf, some_inf) != sub(some_inf, some_inf));
  169. svd.compute(MatrixType::Constant(10,10,some_inf), ComputeFullU | ComputeFullV);
  170. Scalar some_nan = zero<Scalar> () / zero<Scalar> ();
  171. VERIFY(some_nan != some_nan);
  172. svd.compute(MatrixType::Constant(10,10,some_nan), ComputeFullU | ComputeFullV);
  173. MatrixType m = MatrixType::Zero(10,10);
  174. m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_inf;
  175. svd.compute(m, ComputeFullU | ComputeFullV);
  176. m = MatrixType::Zero(10,10);
  177. m(internal::random<int>(0,9), internal::random<int>(0,9)) = some_nan;
  178. svd.compute(m, ComputeFullU | ComputeFullV);
  179. }
  180. template<typename SVD>
  181. void svd_preallocate()
  182. {
  183. Vector3f v(3.f, 2.f, 1.f);
  184. MatrixXf m = v.asDiagonal();
  185. internal::set_is_malloc_allowed(false);
  186. VERIFY_RAISES_ASSERT(VectorXf v(10);)
  187. SVD svd;
  188. internal::set_is_malloc_allowed(true);
  189. svd.compute(m);
  190. VERIFY_IS_APPROX(svd.singularValues(), v);
  191. SVD svd2(3,3);
  192. internal::set_is_malloc_allowed(false);
  193. svd2.compute(m);
  194. internal::set_is_malloc_allowed(true);
  195. VERIFY_IS_APPROX(svd2.singularValues(), v);
  196. VERIFY_RAISES_ASSERT(svd2.matrixU());
  197. VERIFY_RAISES_ASSERT(svd2.matrixV());
  198. svd2.compute(m, ComputeFullU | ComputeFullV);
  199. VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
  200. VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
  201. internal::set_is_malloc_allowed(false);
  202. svd2.compute(m);
  203. internal::set_is_malloc_allowed(true);
  204. SVD svd3(3,3,ComputeFullU|ComputeFullV);
  205. internal::set_is_malloc_allowed(false);
  206. svd2.compute(m);
  207. internal::set_is_malloc_allowed(true);
  208. VERIFY_IS_APPROX(svd2.singularValues(), v);
  209. VERIFY_IS_APPROX(svd2.matrixU(), Matrix3f::Identity());
  210. VERIFY_IS_APPROX(svd2.matrixV(), Matrix3f::Identity());
  211. internal::set_is_malloc_allowed(false);
  212. svd2.compute(m, ComputeFullU|ComputeFullV);
  213. internal::set_is_malloc_allowed(true);
  214. }