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 // 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 template<typename ExpressionType, typename Scalar>
00017 inline void stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
00018 {
00019   Scalar max = bl.cwiseAbs().maxCoeff();
00020   if (max>scale)
00021   {
00022     ssq = ssq * abs2(scale/max);
00023     scale = max;
00024     invScale = Scalar(1)/scale;
00025   }
00026   // TODO if the max is much much smaller than the current scale,
00027   // then we can neglect this sub vector
00028   ssq += (bl*invScale).squaredNorm();
00029 }
00030 }
00031 
00042 template<typename Derived>
00043 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00044 MatrixBase<Derived>::stableNorm() const
00045 {
00046   using std::min;
00047   const Index blockSize = 4096;
00048   RealScalar scale(0);
00049   RealScalar invScale(1);
00050   RealScalar ssq(0); // sum of square
00051   enum {
00052     Alignment = (int(Flags)&DirectAccessBit) || (int(Flags)&AlignedBit) ? 1 : 0
00053   };
00054   Index n = size();
00055   Index bi = internal::first_aligned(derived());
00056   if (bi>0)
00057     internal::stable_norm_kernel(this->head(bi), ssq, scale, invScale);
00058   for (; bi<n; bi+=blockSize)
00059     internal::stable_norm_kernel(this->segment(bi,(min)(blockSize, n - bi)).template forceAlignedAccessIf<Alignment>(), ssq, scale, invScale);
00060   return scale * internal::sqrt(ssq);
00061 }
00062 
00072 template<typename Derived>
00073 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00074 MatrixBase<Derived>::blueNorm() const
00075 {
00076   using std::pow;
00077   using std::min;
00078   using std::max;
00079   static Index nmax = -1;
00080   static RealScalar b1, b2, s1m, s2m, overfl, rbig, relerr;
00081   if(nmax <= 0)
00082   {
00083     int nbig, ibeta, it, iemin, iemax, iexp;
00084     RealScalar abig, eps;
00085     // This program calculates the machine-dependent constants
00086     // bl, b2, slm, s2m, relerr overfl, nmax
00087     // from the "basic" machine-dependent numbers
00088     // nbig, ibeta, it, iemin, iemax, rbig.
00089     // The following define the basic machine-dependent constants.
00090     // For portability, the PORT subprograms "ilmaeh" and "rlmach"
00091     // are used. For any specific computer, each of the assignment
00092     // statements can be replaced
00093     nbig  = (std::numeric_limits<Index>::max)();            // largest integer
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     if (RealScalar(nbig)>abig)  nmax = int(abig);  // largest safe n
00115     else                        nmax = nbig;
00116   }
00117   Index n = size();
00118   RealScalar ab2 = b2 / RealScalar(n);
00119   RealScalar asml = RealScalar(0);
00120   RealScalar amed = RealScalar(0);
00121   RealScalar abig = RealScalar(0);
00122   for(Index j=0; j<n; ++j)
00123   {
00124     RealScalar ax = internal::abs(coeff(j));
00125     if(ax > ab2)     abig += internal::abs2(ax*s2m);
00126     else if(ax < b1) asml += internal::abs2(ax*s1m);
00127     else             amed += internal::abs2(ax);
00128   }
00129   if(abig > RealScalar(0))
00130   {
00131     abig = internal::sqrt(abig);
00132     if(abig > overfl)
00133     {
00134       return rbig;
00135     }
00136     if(amed > RealScalar(0))
00137     {
00138       abig = abig/s2m;
00139       amed = internal::sqrt(amed);
00140     }
00141     else
00142       return abig/s2m;
00143   }
00144   else if(asml > RealScalar(0))
00145   {
00146     if (amed > RealScalar(0))
00147     {
00148       abig = internal::sqrt(amed);
00149       amed = internal::sqrt(asml) / s1m;
00150     }
00151     else
00152       return internal::sqrt(asml)/s1m;
00153   }
00154   else
00155     return internal::sqrt(amed);
00156   asml = (min)(abig, amed);
00157   abig = (max)(abig, amed);
00158   if(asml <= abig*relerr)
00159     return abig;
00160   else
00161     return abig * internal::sqrt(RealScalar(1) + internal::abs2(asml/abig));
00162 }
00163 
00169 template<typename Derived>
00170 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00171 MatrixBase<Derived>::hypotNorm() const
00172 {
00173   return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
00174 }
00175 
00176 } // end namespace Eigen
00177 
00178 #endif // EIGEN_STABLENORM_H


win_eigen
Author(s): Daniel Stonier
autogenerated on Mon Oct 6 2014 12:26:27