1 // This file is part of Eigen, a lightweight C++ template library 2 // for linear algebra. 3 // 4 // Copyright (C) 2009-2014 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> EIGEN_DONT_INLINE T copy(const T& x) 13 { 14 return x; 15 } 16 17 template<typename MatrixType> void stable_norm(const MatrixType& m) 18 { 19 /* this test covers the following files: 20 StableNorm.h 21 */ 22 using std::sqrt; 23 using std::abs; 24 typedef typename MatrixType::Index Index; 25 typedef typename MatrixType::Scalar Scalar; 26 typedef typename NumTraits<Scalar>::Real RealScalar; 27 28 bool complex_real_product_ok = true; 29 30 // Check the basic machine-dependent constants. 31 { 32 int ibeta, it, iemin, iemax; 33 34 ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers 35 it = std::numeric_limits<RealScalar>::digits; // number of base-beta digits in mantissa 36 iemin = std::numeric_limits<RealScalar>::min_exponent; // minimum exponent 37 iemax = std::numeric_limits<RealScalar>::max_exponent; // maximum exponent 38 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 Scalar inf = std::numeric_limits<RealScalar>::infinity(); 43 if(NumTraits<Scalar>::IsComplex && (numext::isnan)(inf*RealScalar(1)) ) 44 { 45 complex_real_product_ok = false; 46 static bool first = true; 47 if(first) 48 std::cerr << "WARNING: compiler mess up complex*real product, " << inf << " * " << 1.0 << " = " << inf*RealScalar(1) << std::endl; 49 first = false; 50 } 51 } 52 53 54 Index rows = m.rows(); 55 Index cols = m.cols(); 56 57 // get a non-zero random factor 58 Scalar factor = internal::random<Scalar>(); 59 while(numext::abs2(factor)<RealScalar(1e-4)) 60 factor = internal::random<Scalar>(); 61 Scalar big = factor * ((std::numeric_limits<RealScalar>::max)() * RealScalar(1e-4)); 62 63 factor = internal::random<Scalar>(); 64 while(numext::abs2(factor)<RealScalar(1e-4)) 65 factor = internal::random<Scalar>(); 66 Scalar small = factor * ((std::numeric_limits<RealScalar>::min)() * RealScalar(1e4)); 67 68 MatrixType vzero = MatrixType::Zero(rows, cols), 69 vrand = MatrixType::Random(rows, cols), 70 vbig(rows, cols), 71 vsmall(rows,cols); 72 73 vbig.fill(big); 74 vsmall.fill(small); 75 76 VERIFY_IS_MUCH_SMALLER_THAN(vzero.norm(), static_cast<RealScalar>(1)); 77 VERIFY_IS_APPROX(vrand.stableNorm(), vrand.norm()); 78 VERIFY_IS_APPROX(vrand.blueNorm(), vrand.norm()); 79 VERIFY_IS_APPROX(vrand.hypotNorm(), vrand.norm()); 80 81 RealScalar size = static_cast<RealScalar>(m.size()); 82 83 // test numext::isfinite 84 VERIFY(!(numext::isfinite)( std::numeric_limits<RealScalar>::infinity())); 85 VERIFY(!(numext::isfinite)(sqrt(-abs(big)))); 86 87 // test overflow 88 VERIFY((numext::isfinite)(sqrt(size)*abs(big))); 89 VERIFY_IS_NOT_APPROX(sqrt(copy(vbig.squaredNorm())), abs(sqrt(size)*big)); // here the default norm must fail 90 VERIFY_IS_APPROX(vbig.stableNorm(), sqrt(size)*abs(big)); 91 VERIFY_IS_APPROX(vbig.blueNorm(), sqrt(size)*abs(big)); 92 VERIFY_IS_APPROX(vbig.hypotNorm(), sqrt(size)*abs(big)); 93 94 // test underflow 95 VERIFY((numext::isfinite)(sqrt(size)*abs(small))); 96 VERIFY_IS_NOT_APPROX(sqrt(copy(vsmall.squaredNorm())), abs(sqrt(size)*small)); // here the default norm must fail 97 VERIFY_IS_APPROX(vsmall.stableNorm(), sqrt(size)*abs(small)); 98 VERIFY_IS_APPROX(vsmall.blueNorm(), sqrt(size)*abs(small)); 99 VERIFY_IS_APPROX(vsmall.hypotNorm(), sqrt(size)*abs(small)); 100 101 // Test compilation of cwise() version 102 VERIFY_IS_APPROX(vrand.colwise().stableNorm(), vrand.colwise().norm()); 103 VERIFY_IS_APPROX(vrand.colwise().blueNorm(), vrand.colwise().norm()); 104 VERIFY_IS_APPROX(vrand.colwise().hypotNorm(), vrand.colwise().norm()); 105 VERIFY_IS_APPROX(vrand.rowwise().stableNorm(), vrand.rowwise().norm()); 106 VERIFY_IS_APPROX(vrand.rowwise().blueNorm(), vrand.rowwise().norm()); 107 VERIFY_IS_APPROX(vrand.rowwise().hypotNorm(), vrand.rowwise().norm()); 108 109 // test NaN, +inf, -inf 110 MatrixType v; 111 Index i = internal::random<Index>(0,rows-1); 112 Index j = internal::random<Index>(0,cols-1); 113 114 // NaN 115 { 116 v = vrand; 117 v(i,j) = std::numeric_limits<RealScalar>::quiet_NaN(); 118 VERIFY(!(numext::isfinite)(v.squaredNorm())); VERIFY((numext::isnan)(v.squaredNorm())); 119 VERIFY(!(numext::isfinite)(v.norm())); VERIFY((numext::isnan)(v.norm())); 120 VERIFY(!(numext::isfinite)(v.stableNorm())); VERIFY((numext::isnan)(v.stableNorm())); 121 VERIFY(!(numext::isfinite)(v.blueNorm())); VERIFY((numext::isnan)(v.blueNorm())); 122 VERIFY(!(numext::isfinite)(v.hypotNorm())); VERIFY((numext::isnan)(v.hypotNorm())); 123 } 124 125 // +inf 126 { 127 v = vrand; 128 v(i,j) = std::numeric_limits<RealScalar>::infinity(); 129 VERIFY(!(numext::isfinite)(v.squaredNorm())); VERIFY(isPlusInf(v.squaredNorm())); 130 VERIFY(!(numext::isfinite)(v.norm())); VERIFY(isPlusInf(v.norm())); 131 VERIFY(!(numext::isfinite)(v.stableNorm())); 132 if(complex_real_product_ok){ 133 VERIFY(isPlusInf(v.stableNorm())); 134 } 135 VERIFY(!(numext::isfinite)(v.blueNorm())); VERIFY(isPlusInf(v.blueNorm())); 136 VERIFY(!(numext::isfinite)(v.hypotNorm())); VERIFY(isPlusInf(v.hypotNorm())); 137 } 138 139 // -inf 140 { 141 v = vrand; 142 v(i,j) = -std::numeric_limits<RealScalar>::infinity(); 143 VERIFY(!(numext::isfinite)(v.squaredNorm())); VERIFY(isPlusInf(v.squaredNorm())); 144 VERIFY(!(numext::isfinite)(v.norm())); VERIFY(isPlusInf(v.norm())); 145 VERIFY(!(numext::isfinite)(v.stableNorm())); 146 if(complex_real_product_ok) { 147 VERIFY(isPlusInf(v.stableNorm())); 148 } 149 VERIFY(!(numext::isfinite)(v.blueNorm())); VERIFY(isPlusInf(v.blueNorm())); 150 VERIFY(!(numext::isfinite)(v.hypotNorm())); VERIFY(isPlusInf(v.hypotNorm())); 151 } 152 153 // mix 154 { 155 Index i2 = internal::random<Index>(0,rows-1); 156 Index j2 = internal::random<Index>(0,cols-1); 157 v = vrand; 158 v(i,j) = -std::numeric_limits<RealScalar>::infinity(); 159 v(i2,j2) = std::numeric_limits<RealScalar>::quiet_NaN(); 160 VERIFY(!(numext::isfinite)(v.squaredNorm())); VERIFY((numext::isnan)(v.squaredNorm())); 161 VERIFY(!(numext::isfinite)(v.norm())); VERIFY((numext::isnan)(v.norm())); 162 VERIFY(!(numext::isfinite)(v.stableNorm())); VERIFY((numext::isnan)(v.stableNorm())); 163 VERIFY(!(numext::isfinite)(v.blueNorm())); VERIFY((numext::isnan)(v.blueNorm())); 164 VERIFY(!(numext::isfinite)(v.hypotNorm())); VERIFY((numext::isnan)(v.hypotNorm())); 165 } 166 167 // stableNormalize[d] 168 { 169 VERIFY_IS_APPROX(vrand.stableNormalized(), vrand.normalized()); 170 MatrixType vcopy(vrand); 171 vcopy.stableNormalize(); 172 VERIFY_IS_APPROX(vcopy, vrand.normalized()); 173 VERIFY_IS_APPROX((vrand.stableNormalized()).norm(), RealScalar(1)); 174 VERIFY_IS_APPROX(vcopy.norm(), RealScalar(1)); 175 VERIFY_IS_APPROX((vbig.stableNormalized()).norm(), RealScalar(1)); 176 VERIFY_IS_APPROX((vsmall.stableNormalized()).norm(), RealScalar(1)); 177 RealScalar big_scaling = ((std::numeric_limits<RealScalar>::max)() * RealScalar(1e-4)); 178 VERIFY_IS_APPROX(vbig/big_scaling, (vbig.stableNorm() * vbig.stableNormalized()).eval()/big_scaling); 179 VERIFY_IS_APPROX(vsmall, vsmall.stableNorm() * vsmall.stableNormalized()); 180 } 181 } 182 183 void test_stable_norm() 184 { 185 for(int i = 0; i < g_repeat; i++) { 186 CALL_SUBTEST_1( stable_norm(Matrix<float, 1, 1>()) ); 187 CALL_SUBTEST_2( stable_norm(Vector4d()) ); 188 CALL_SUBTEST_3( stable_norm(VectorXd(internal::random<int>(10,2000))) ); 189 CALL_SUBTEST_4( stable_norm(VectorXf(internal::random<int>(10,2000))) ); 190 CALL_SUBTEST_5( stable_norm(VectorXcd(internal::random<int>(10,2000))) ); 191 } 192 } 193