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.

110 lines
4.3 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "main.h"
  10. template<typename T> bool isNotNaN(const T& x)
  11. {
  12. return x==x;
  13. }
  14. // workaround aggressive optimization in ICC
  15. template<typename T> EIGEN_DONT_INLINE T sub(T a, T b) { return a - b; }
  16. template<typename T> bool isFinite(const T& x)
  17. {
  18. return isNotNaN(sub(x,x));
  19. }
  20. template<typename T> EIGEN_DONT_INLINE T copy(const T& x)
  21. {
  22. return x;
  23. }
  24. template<typename MatrixType> void stable_norm(const MatrixType& m)
  25. {
  26. /* this test covers the following files:
  27. StableNorm.h
  28. */
  29. typedef typename MatrixType::Index Index;
  30. typedef typename MatrixType::Scalar Scalar;
  31. typedef typename NumTraits<Scalar>::Real RealScalar;
  32. // Check the basic machine-dependent constants.
  33. {
  34. int ibeta, it, iemin, iemax;
  35. ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers
  36. it = std::numeric_limits<RealScalar>::digits; // number of base-beta digits in mantissa
  37. iemin = std::numeric_limits<RealScalar>::min_exponent; // minimum exponent
  38. iemax = std::numeric_limits<RealScalar>::max_exponent; // maximum exponent
  39. VERIFY( (!(iemin > 1 - 2*it || 1+it>iemax || (it==2 && ibeta<5) || (it<=4 && ibeta <= 3 ) || it<2))
  40. && "the stable norm algorithm cannot be guaranteed on this computer");
  41. }
  42. Index rows = m.rows();
  43. Index cols = m.cols();
  44. Scalar big = internal::random<Scalar>() * ((std::numeric_limits<RealScalar>::max)() * RealScalar(1e-4));
  45. Scalar small = internal::random<Scalar>() * ((std::numeric_limits<RealScalar>::min)() * RealScalar(1e4));
  46. MatrixType vzero = MatrixType::Zero(rows, cols),
  47. vrand = MatrixType::Random(rows, cols),
  48. vbig(rows, cols),
  49. vsmall(rows,cols);
  50. vbig.fill(big);
  51. vsmall.fill(small);
  52. VERIFY_IS_MUCH_SMALLER_THAN(vzero.norm(), static_cast<RealScalar>(1));
  53. VERIFY_IS_APPROX(vrand.stableNorm(), vrand.norm());
  54. VERIFY_IS_APPROX(vrand.blueNorm(), vrand.norm());
  55. VERIFY_IS_APPROX(vrand.hypotNorm(), vrand.norm());
  56. RealScalar size = static_cast<RealScalar>(m.size());
  57. // test isFinite
  58. VERIFY(!isFinite( std::numeric_limits<RealScalar>::infinity()));
  59. VERIFY(!isFinite(internal::sqrt(-internal::abs(big))));
  60. // test overflow
  61. VERIFY(isFinite(internal::sqrt(size)*internal::abs(big)));
  62. VERIFY_IS_NOT_APPROX(internal::sqrt(copy(vbig.squaredNorm())), internal::abs(internal::sqrt(size)*big)); // here the default norm must fail
  63. VERIFY_IS_APPROX(vbig.stableNorm(), internal::sqrt(size)*internal::abs(big));
  64. VERIFY_IS_APPROX(vbig.blueNorm(), internal::sqrt(size)*internal::abs(big));
  65. VERIFY_IS_APPROX(vbig.hypotNorm(), internal::sqrt(size)*internal::abs(big));
  66. // test underflow
  67. VERIFY(isFinite(internal::sqrt(size)*internal::abs(small)));
  68. VERIFY_IS_NOT_APPROX(internal::sqrt(copy(vsmall.squaredNorm())), internal::abs(internal::sqrt(size)*small)); // here the default norm must fail
  69. VERIFY_IS_APPROX(vsmall.stableNorm(), internal::sqrt(size)*internal::abs(small));
  70. VERIFY_IS_APPROX(vsmall.blueNorm(), internal::sqrt(size)*internal::abs(small));
  71. VERIFY_IS_APPROX(vsmall.hypotNorm(), internal::sqrt(size)*internal::abs(small));
  72. // Test compilation of cwise() version
  73. VERIFY_IS_APPROX(vrand.colwise().stableNorm(), vrand.colwise().norm());
  74. VERIFY_IS_APPROX(vrand.colwise().blueNorm(), vrand.colwise().norm());
  75. VERIFY_IS_APPROX(vrand.colwise().hypotNorm(), vrand.colwise().norm());
  76. VERIFY_IS_APPROX(vrand.rowwise().stableNorm(), vrand.rowwise().norm());
  77. VERIFY_IS_APPROX(vrand.rowwise().blueNorm(), vrand.rowwise().norm());
  78. VERIFY_IS_APPROX(vrand.rowwise().hypotNorm(), vrand.rowwise().norm());
  79. }
  80. void test_stable_norm()
  81. {
  82. for(int i = 0; i < g_repeat; i++) {
  83. CALL_SUBTEST_1( stable_norm(Matrix<float, 1, 1>()) );
  84. CALL_SUBTEST_2( stable_norm(Vector4d()) );
  85. CALL_SUBTEST_3( stable_norm(VectorXd(internal::random<int>(10,2000))) );
  86. CALL_SUBTEST_4( stable_norm(VectorXf(internal::random<int>(10,2000))) );
  87. CALL_SUBTEST_5( stable_norm(VectorXcd(internal::random<int>(10,2000))) );
  88. }
  89. }