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.4 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. using std::sqrt;
  30. using std::abs;
  31. typedef typename MatrixType::Index Index;
  32. typedef typename MatrixType::Scalar Scalar;
  33. typedef typename NumTraits<Scalar>::Real RealScalar;
  34. // Check the basic machine-dependent constants.
  35. {
  36. int ibeta, it, iemin, iemax;
  37. ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers
  38. it = std::numeric_limits<RealScalar>::digits; // number of base-beta digits in mantissa
  39. iemin = std::numeric_limits<RealScalar>::min_exponent; // minimum exponent
  40. iemax = std::numeric_limits<RealScalar>::max_exponent; // maximum exponent
  41. VERIFY( (!(iemin > 1 - 2*it || 1+it>iemax || (it==2 && ibeta<5) || (it<=4 && ibeta <= 3 ) || it<2))
  42. && "the stable norm algorithm cannot be guaranteed on this computer");
  43. }
  44. Index rows = m.rows();
  45. Index cols = m.cols();
  46. // get a non-zero random factor
  47. Scalar factor = internal::random<Scalar>();
  48. while(numext::abs2(factor)<RealScalar(1e-4))
  49. factor = internal::random<Scalar>();
  50. Scalar big = factor * ((std::numeric_limits<RealScalar>::max)() * RealScalar(1e-4));
  51. factor = internal::random<Scalar>();
  52. while(numext::abs2(factor)<RealScalar(1e-4))
  53. factor = internal::random<Scalar>();
  54. Scalar small = factor * ((std::numeric_limits<RealScalar>::min)() * RealScalar(1e4));
  55. MatrixType vzero = MatrixType::Zero(rows, cols),
  56. vrand = MatrixType::Random(rows, cols),
  57. vbig(rows, cols),
  58. vsmall(rows,cols);
  59. vbig.fill(big);
  60. vsmall.fill(small);
  61. VERIFY_IS_MUCH_SMALLER_THAN(vzero.norm(), static_cast<RealScalar>(1));
  62. VERIFY_IS_APPROX(vrand.stableNorm(), vrand.norm());
  63. VERIFY_IS_APPROX(vrand.blueNorm(), vrand.norm());
  64. VERIFY_IS_APPROX(vrand.hypotNorm(), vrand.norm());
  65. RealScalar size = static_cast<RealScalar>(m.size());
  66. // test isFinite
  67. VERIFY(!isFinite( std::numeric_limits<RealScalar>::infinity()));
  68. VERIFY(!isFinite(sqrt(-abs(big))));
  69. // test overflow
  70. VERIFY(isFinite(sqrt(size)*abs(big)));
  71. VERIFY_IS_NOT_APPROX(sqrt(copy(vbig.squaredNorm())), abs(sqrt(size)*big)); // here the default norm must fail
  72. VERIFY_IS_APPROX(vbig.stableNorm(), sqrt(size)*abs(big));
  73. VERIFY_IS_APPROX(vbig.blueNorm(), sqrt(size)*abs(big));
  74. VERIFY_IS_APPROX(vbig.hypotNorm(), sqrt(size)*abs(big));
  75. // test underflow
  76. VERIFY(isFinite(sqrt(size)*abs(small)));
  77. VERIFY_IS_NOT_APPROX(sqrt(copy(vsmall.squaredNorm())), abs(sqrt(size)*small)); // here the default norm must fail
  78. VERIFY_IS_APPROX(vsmall.stableNorm(), sqrt(size)*abs(small));
  79. VERIFY_IS_APPROX(vsmall.blueNorm(), sqrt(size)*abs(small));
  80. VERIFY_IS_APPROX(vsmall.hypotNorm(), sqrt(size)*abs(small));
  81. // Test compilation of cwise() version
  82. VERIFY_IS_APPROX(vrand.colwise().stableNorm(), vrand.colwise().norm());
  83. VERIFY_IS_APPROX(vrand.colwise().blueNorm(), vrand.colwise().norm());
  84. VERIFY_IS_APPROX(vrand.colwise().hypotNorm(), vrand.colwise().norm());
  85. VERIFY_IS_APPROX(vrand.rowwise().stableNorm(), vrand.rowwise().norm());
  86. VERIFY_IS_APPROX(vrand.rowwise().blueNorm(), vrand.rowwise().norm());
  87. VERIFY_IS_APPROX(vrand.rowwise().hypotNorm(), vrand.rowwise().norm());
  88. }
  89. void test_stable_norm()
  90. {
  91. for(int i = 0; i < g_repeat; i++) {
  92. CALL_SUBTEST_1( stable_norm(Matrix<float, 1, 1>()) );
  93. CALL_SUBTEST_2( stable_norm(Vector4d()) );
  94. CALL_SUBTEST_3( stable_norm(VectorXd(internal::random<int>(10,2000))) );
  95. CALL_SUBTEST_4( stable_norm(VectorXf(internal::random<int>(10,2000))) );
  96. CALL_SUBTEST_5( stable_norm(VectorXcd(internal::random<int>(10,2000))) );
  97. }
  98. }