StableNorm.h
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 //
00006 // This Source Code Form is subject to the terms of the Mozilla
00007 // Public License v. 2.0. If a copy of the MPL was not distributed
00008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
00009 
00010 #ifndef EIGEN_STABLENORM_H
00011 #define EIGEN_STABLENORM_H
00012 
00013 namespace Eigen { 
00014 
00015 namespace internal {
00016 
00017 template<typename ExpressionType, typename Scalar>
00018 inline void stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
00019 {
00020   Scalar max = bl.cwiseAbs().maxCoeff();
00021   if (max>scale)
00022   {
00023     ssq = ssq * abs2(scale/max);
00024     scale = max;
00025     invScale = Scalar(1)/scale;
00026   }
00027   // TODO if the max is much much smaller than the current scale,
00028   // then we can neglect this sub vector
00029   ssq += (bl*invScale).squaredNorm();
00030 }
00031 }
00032 
00043 template<typename Derived>
00044 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00045 MatrixBase<Derived>::stableNorm() const
00046 {
00047   using std::min;
00048   const Index blockSize = 4096;
00049   RealScalar scale(0);
00050   RealScalar invScale(1);
00051   RealScalar ssq(0); // sum of square
00052   enum {
00053     Alignment = (int(Flags)&DirectAccessBit) || (int(Flags)&AlignedBit) ? 1 : 0
00054   };
00055   Index n = size();
00056   Index bi = internal::first_aligned(derived());
00057   if (bi>0)
00058     internal::stable_norm_kernel(this->head(bi), ssq, scale, invScale);
00059   for (; bi<n; bi+=blockSize)
00060     internal::stable_norm_kernel(this->segment(bi,(min)(blockSize, n - bi)).template forceAlignedAccessIf<Alignment>(), ssq, scale, invScale);
00061   return scale * internal::sqrt(ssq);
00062 }
00063 
00073 template<typename Derived>
00074 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00075 MatrixBase<Derived>::blueNorm() const
00076 {
00077   using std::pow;
00078   using std::min;
00079   using std::max;
00080   static bool initialized = false;
00081   static RealScalar b1, b2, s1m, s2m, overfl, rbig, relerr;
00082   if(!initialized)
00083   {
00084     int ibeta, it, iemin, iemax, iexp;
00085     RealScalar abig, eps;
00086     // This program calculates the machine-dependent constants
00087     // bl, b2, slm, s2m, relerr overfl
00088     // from the "basic" machine-dependent numbers
00089     // ibeta, it, iemin, iemax, rbig.
00090     // The following define the basic machine-dependent constants.
00091     // For portability, the PORT subprograms "ilmaeh" and "rlmach"
00092     // are used. For any specific computer, each of the assignment
00093     // statements can be replaced
00094     ibeta = std::numeric_limits<RealScalar>::radix;         // base for floating-point numbers
00095     it    = std::numeric_limits<RealScalar>::digits;        // number of base-beta digits in mantissa
00096     iemin = std::numeric_limits<RealScalar>::min_exponent;  // minimum exponent
00097     iemax = std::numeric_limits<RealScalar>::max_exponent;  // maximum exponent
00098     rbig  = (std::numeric_limits<RealScalar>::max)();         // largest floating-point number
00099 
00100     iexp  = -((1-iemin)/2);
00101     b1    = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));  // lower boundary of midrange
00102     iexp  = (iemax + 1 - it)/2;
00103     b2    = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));   // upper boundary of midrange
00104 
00105     iexp  = (2-iemin)/2;
00106     s1m   = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));   // scaling factor for lower range
00107     iexp  = - ((iemax+it)/2);
00108     s2m   = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));   // scaling factor for upper range
00109 
00110     overfl  = rbig*s2m;             // overflow boundary for abig
00111     eps     = RealScalar(pow(double(ibeta), 1-it));
00112     relerr  = internal::sqrt(eps);         // tolerance for neglecting asml
00113     abig    = RealScalar(1.0/eps - 1.0);
00114     initialized = true;
00115   }
00116   Index n = size();
00117   RealScalar ab2 = b2 / RealScalar(n);
00118   RealScalar asml = RealScalar(0);
00119   RealScalar amed = RealScalar(0);
00120   RealScalar abig = RealScalar(0);
00121   for(Index j=0; j<n; ++j)
00122   {
00123     RealScalar ax = internal::abs(coeff(j));
00124     if(ax > ab2)     abig += internal::abs2(ax*s2m);
00125     else if(ax < b1) asml += internal::abs2(ax*s1m);
00126     else             amed += internal::abs2(ax);
00127   }
00128   if(abig > RealScalar(0))
00129   {
00130     abig = internal::sqrt(abig);
00131     if(abig > overfl)
00132     {
00133       return rbig;
00134     }
00135     if(amed > RealScalar(0))
00136     {
00137       abig = abig/s2m;
00138       amed = internal::sqrt(amed);
00139     }
00140     else
00141       return abig/s2m;
00142   }
00143   else if(asml > RealScalar(0))
00144   {
00145     if (amed > RealScalar(0))
00146     {
00147       abig = internal::sqrt(amed);
00148       amed = internal::sqrt(asml) / s1m;
00149     }
00150     else
00151       return internal::sqrt(asml)/s1m;
00152   }
00153   else
00154     return internal::sqrt(amed);
00155   asml = (min)(abig, amed);
00156   abig = (max)(abig, amed);
00157   if(asml <= abig*relerr)
00158     return abig;
00159   else
00160     return abig * internal::sqrt(RealScalar(1) + internal::abs2(asml/abig));
00161 }
00162 
00168 template<typename Derived>
00169 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00170 MatrixBase<Derived>::hypotNorm() const
00171 {
00172   return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
00173 }
00174 
00175 } // end namespace Eigen
00176 
00177 #endif // EIGEN_STABLENORM_H