Home | History | Annotate | Download | only in test
      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 (at) 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 
     10 #include "main.h"
     11 
     12 template<typename T> bool isNotNaN(const T& x)
     13 {
     14   return x==x;
     15 }
     16 
     17 // workaround aggressive optimization in ICC
     18 template<typename T> EIGEN_DONT_INLINE  T sub(T a, T b) { return a - b; }
     19 
     20 template<typename T> bool isFinite(const T& x)
     21 {
     22   return isNotNaN(sub(x,x));
     23 }
     24 
     25 template<typename T> EIGEN_DONT_INLINE T copy(const T& x)
     26 {
     27   return x;
     28 }
     29 
     30 template<typename MatrixType> void stable_norm(const MatrixType& m)
     31 {
     32   /* this test covers the following files:
     33      StableNorm.h
     34   */
     35   typedef typename MatrixType::Index Index;
     36   typedef typename MatrixType::Scalar Scalar;
     37   typedef typename NumTraits<Scalar>::Real RealScalar;
     38 
     39   // Check the basic machine-dependent constants.
     40   {
     41     int ibeta, it, iemin, iemax;
     42 
     43     ibeta = std::numeric_limits<RealScalar>::radix;         // base for floating-point numbers
     44     it    = std::numeric_limits<RealScalar>::digits;        // number of base-beta digits in mantissa
     45     iemin = std::numeric_limits<RealScalar>::min_exponent;  // minimum exponent
     46     iemax = std::numeric_limits<RealScalar>::max_exponent;  // maximum exponent
     47 
     48     VERIFY( (!(iemin > 1 - 2*it || 1+it>iemax || (it==2 && ibeta<5) || (it<=4 && ibeta <= 3 ) || it<2))
     49            && "the stable norm algorithm cannot be guaranteed on this computer");
     50   }
     51 
     52 
     53   Index rows = m.rows();
     54   Index cols = m.cols();
     55 
     56   Scalar big = internal::random<Scalar>() * ((std::numeric_limits<RealScalar>::max)() * RealScalar(1e-4));
     57   Scalar small = internal::random<Scalar>() * ((std::numeric_limits<RealScalar>::min)() * RealScalar(1e4));
     58 
     59   MatrixType  vzero = MatrixType::Zero(rows, cols),
     60               vrand = MatrixType::Random(rows, cols),
     61               vbig(rows, cols),
     62               vsmall(rows,cols);
     63 
     64   vbig.fill(big);
     65   vsmall.fill(small);
     66 
     67   VERIFY_IS_MUCH_SMALLER_THAN(vzero.norm(), static_cast<RealScalar>(1));
     68   VERIFY_IS_APPROX(vrand.stableNorm(),      vrand.norm());
     69   VERIFY_IS_APPROX(vrand.blueNorm(),        vrand.norm());
     70   VERIFY_IS_APPROX(vrand.hypotNorm(),       vrand.norm());
     71 
     72   RealScalar size = static_cast<RealScalar>(m.size());
     73 
     74   // test isFinite
     75   VERIFY(!isFinite( std::numeric_limits<RealScalar>::infinity()));
     76   VERIFY(!isFinite(internal::sqrt(-internal::abs(big))));
     77 
     78   // test overflow
     79   VERIFY(isFinite(internal::sqrt(size)*internal::abs(big)));
     80   VERIFY_IS_NOT_APPROX(internal::sqrt(copy(vbig.squaredNorm())),   internal::abs(internal::sqrt(size)*big)); // here the default norm must fail
     81   VERIFY_IS_APPROX(vbig.stableNorm(), internal::sqrt(size)*internal::abs(big));
     82   VERIFY_IS_APPROX(vbig.blueNorm(),   internal::sqrt(size)*internal::abs(big));
     83   VERIFY_IS_APPROX(vbig.hypotNorm(),  internal::sqrt(size)*internal::abs(big));
     84 
     85   // test underflow
     86   VERIFY(isFinite(internal::sqrt(size)*internal::abs(small)));
     87   VERIFY_IS_NOT_APPROX(internal::sqrt(copy(vsmall.squaredNorm())),   internal::abs(internal::sqrt(size)*small)); // here the default norm must fail
     88   VERIFY_IS_APPROX(vsmall.stableNorm(), internal::sqrt(size)*internal::abs(small));
     89   VERIFY_IS_APPROX(vsmall.blueNorm(),   internal::sqrt(size)*internal::abs(small));
     90   VERIFY_IS_APPROX(vsmall.hypotNorm(),  internal::sqrt(size)*internal::abs(small));
     91 
     92 // Test compilation of cwise() version
     93   VERIFY_IS_APPROX(vrand.colwise().stableNorm(),      vrand.colwise().norm());
     94   VERIFY_IS_APPROX(vrand.colwise().blueNorm(),        vrand.colwise().norm());
     95   VERIFY_IS_APPROX(vrand.colwise().hypotNorm(),       vrand.colwise().norm());
     96   VERIFY_IS_APPROX(vrand.rowwise().stableNorm(),      vrand.rowwise().norm());
     97   VERIFY_IS_APPROX(vrand.rowwise().blueNorm(),        vrand.rowwise().norm());
     98   VERIFY_IS_APPROX(vrand.rowwise().hypotNorm(),       vrand.rowwise().norm());
     99 }
    100 
    101 void test_stable_norm()
    102 {
    103   for(int i = 0; i < g_repeat; i++) {
    104     CALL_SUBTEST_1( stable_norm(Matrix<float, 1, 1>()) );
    105     CALL_SUBTEST_2( stable_norm(Vector4d()) );
    106     CALL_SUBTEST_3( stable_norm(VectorXd(internal::random<int>(10,2000))) );
    107     CALL_SUBTEST_4( stable_norm(VectorXf(internal::random<int>(10,2000))) );
    108     CALL_SUBTEST_5( stable_norm(VectorXcd(internal::random<int>(10,2000))) );
    109   }
    110 }
    111