StableNorm.h
Go to the documentation of this file.
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 // Eigen is free software; you can redistribute it and/or
00007 // modify it under the terms of the GNU Lesser General Public
00008 // License as published by the Free Software Foundation; either
00009 // version 3 of the License, or (at your option) any later version.
00010 //
00011 // Alternatively, you can redistribute it and/or
00012 // modify it under the terms of the GNU General Public License as
00013 // published by the Free Software Foundation; either version 2 of
00014 // the License, or (at your option) any later version.
00015 //
00016 // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
00017 // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
00018 // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
00019 // GNU General Public License for more details.
00020 //
00021 // You should have received a copy of the GNU Lesser General Public
00022 // License and a copy of the GNU General Public License along with
00023 // Eigen. If not, see <http://www.gnu.org/licenses/>.
00024 
00025 #ifndef EIGEN_STABLENORM_H
00026 #define EIGEN_STABLENORM_H
00027 
00028 namespace internal {
00029 template<typename ExpressionType, typename Scalar>
00030 inline void stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
00031 {
00032   Scalar max = bl.cwiseAbs().maxCoeff();
00033   if (max>scale)
00034   {
00035     ssq = ssq * abs2(scale/max);
00036     scale = max;
00037     invScale = Scalar(1)/scale;
00038   }
00039   // TODO if the max is much much smaller than the current scale,
00040   // then we can neglect this sub vector
00041   ssq += (bl*invScale).squaredNorm();
00042 }
00043 }
00044 
00055 template<typename Derived>
00056 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00057 MatrixBase<Derived>::stableNorm() const
00058 {
00059   using std::min;
00060   const Index blockSize = 4096;
00061   RealScalar scale = 0;
00062   RealScalar invScale = 1;
00063   RealScalar ssq = 0; // sum of square
00064   enum {
00065     Alignment = (int(Flags)&DirectAccessBit) || (int(Flags)&AlignedBit) ? 1 : 0
00066   };
00067   Index n = size();
00068   Index bi = internal::first_aligned(derived());
00069   if (bi>0)
00070     internal::stable_norm_kernel(this->head(bi), ssq, scale, invScale);
00071   for (; bi<n; bi+=blockSize)
00072     internal::stable_norm_kernel(this->segment(bi,min(blockSize, n - bi)).template forceAlignedAccessIf<Alignment>(), ssq, scale, invScale);
00073   return scale * internal::sqrt(ssq);
00074 }
00075 
00085 template<typename Derived>
00086 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00087 MatrixBase<Derived>::blueNorm() const
00088 {
00089   using std::pow;
00090   using std::min;
00091   using std::max;
00092   static Index nmax = -1;
00093   static RealScalar b1, b2, s1m, s2m, overfl, rbig, relerr;
00094   if(nmax <= 0)
00095   {
00096     int nbig, ibeta, it, iemin, iemax, iexp;
00097     RealScalar abig, eps;
00098     // This program calculates the machine-dependent constants
00099     // bl, b2, slm, s2m, relerr overfl, nmax
00100     // from the "basic" machine-dependent numbers
00101     // nbig, ibeta, it, iemin, iemax, rbig.
00102     // The following define the basic machine-dependent constants.
00103     // For portability, the PORT subprograms "ilmaeh" and "rlmach"
00104     // are used. For any specific computer, each of the assignment
00105     // statements can be replaced
00106     nbig  = std::numeric_limits<Index>::max();            // largest integer
00107     ibeta = std::numeric_limits<RealScalar>::radix;         // base for floating-point numbers
00108     it    = std::numeric_limits<RealScalar>::digits;        // number of base-beta digits in mantissa
00109     iemin = std::numeric_limits<RealScalar>::min_exponent;  // minimum exponent
00110     iemax = std::numeric_limits<RealScalar>::max_exponent;  // maximum exponent
00111     rbig  = std::numeric_limits<RealScalar>::max();         // largest floating-point number
00112 
00113     iexp  = -((1-iemin)/2);
00114     b1    = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));  // lower boundary of midrange
00115     iexp  = (iemax + 1 - it)/2;
00116     b2    = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));   // upper boundary of midrange
00117 
00118     iexp  = (2-iemin)/2;
00119     s1m   = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));   // scaling factor for lower range
00120     iexp  = - ((iemax+it)/2);
00121     s2m   = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));   // scaling factor for upper range
00122 
00123     overfl  = rbig*s2m;             // overflow boundary for abig
00124     eps     = RealScalar(pow(double(ibeta), 1-it));
00125     relerr  = internal::sqrt(eps);         // tolerance for neglecting asml
00126     abig    = RealScalar(1.0/eps - 1.0);
00127     if (RealScalar(nbig)>abig)  nmax = int(abig);  // largest safe n
00128     else                        nmax = nbig;
00129   }
00130   Index n = size();
00131   RealScalar ab2 = b2 / RealScalar(n);
00132   RealScalar asml = RealScalar(0);
00133   RealScalar amed = RealScalar(0);
00134   RealScalar abig = RealScalar(0);
00135   for(Index j=0; j<n; ++j)
00136   {
00137     RealScalar ax = internal::abs(coeff(j));
00138     if(ax > ab2)     abig += internal::abs2(ax*s2m);
00139     else if(ax < b1) asml += internal::abs2(ax*s1m);
00140     else             amed += internal::abs2(ax);
00141   }
00142   if(abig > RealScalar(0))
00143   {
00144     abig = internal::sqrt(abig);
00145     if(abig > overfl)
00146     {
00147       eigen_assert(false && "overflow");
00148       return rbig;
00149     }
00150     if(amed > RealScalar(0))
00151     {
00152       abig = abig/s2m;
00153       amed = internal::sqrt(amed);
00154     }
00155     else
00156       return abig/s2m;
00157   }
00158   else if(asml > RealScalar(0))
00159   {
00160     if (amed > RealScalar(0))
00161     {
00162       abig = internal::sqrt(amed);
00163       amed = internal::sqrt(asml) / s1m;
00164     }
00165     else
00166       return internal::sqrt(asml)/s1m;
00167   }
00168   else
00169     return internal::sqrt(amed);
00170   asml = min(abig, amed);
00171   abig = max(abig, amed);
00172   if(asml <= abig*relerr)
00173     return abig;
00174   else
00175     return abig * internal::sqrt(RealScalar(1) + internal::abs2(asml/abig));
00176 }
00177 
00183 template<typename Derived>
00184 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00185 MatrixBase<Derived>::hypotNorm() const
00186 {
00187   return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
00188 }
00189 
00190 #endif // EIGEN_STABLENORM_H


re_vision
Author(s): Dorian Galvez-Lopez
autogenerated on Sun Jan 5 2014 11:32:59