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 
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 * numext::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 template<typename Derived>
00033 inline typename NumTraits<typename traits<Derived>::Scalar>::Real
00034 blueNorm_impl(const EigenBase<Derived>& _vec)
00035 {
00036   typedef typename Derived::RealScalar RealScalar;  
00037   typedef typename Derived::Index Index;
00038   using std::pow;
00039   using std::min;
00040   using std::max;
00041   using std::sqrt;
00042   using std::abs;
00043   const Derived& vec(_vec.derived());
00044   static bool initialized = false;
00045   static RealScalar b1, b2, s1m, s2m, overfl, rbig, relerr;
00046   if(!initialized)
00047   {
00048     int ibeta, it, iemin, iemax, iexp;
00049     RealScalar eps;
00050     // This program calculates the machine-dependent constants
00051     // bl, b2, slm, s2m, relerr overfl
00052     // from the "basic" machine-dependent numbers
00053     // nbig, ibeta, it, iemin, iemax, rbig.
00054     // The following define the basic machine-dependent constants.
00055     // For portability, the PORT subprograms "ilmaeh" and "rlmach"
00056     // are used. For any specific computer, each of the assignment
00057     // statements can be replaced
00058     ibeta = std::numeric_limits<RealScalar>::radix;                 // base for floating-point numbers
00059     it    = std::numeric_limits<RealScalar>::digits;                // number of base-beta digits in mantissa
00060     iemin = std::numeric_limits<RealScalar>::min_exponent;          // minimum exponent
00061     iemax = std::numeric_limits<RealScalar>::max_exponent;          // maximum exponent
00062     rbig  = (std::numeric_limits<RealScalar>::max)();               // largest floating-point number
00063 
00064     iexp  = -((1-iemin)/2);
00065     b1    = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));    // lower boundary of midrange
00066     iexp  = (iemax + 1 - it)/2;
00067     b2    = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));    // upper boundary of midrange
00068 
00069     iexp  = (2-iemin)/2;
00070     s1m   = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));    // scaling factor for lower range
00071     iexp  = - ((iemax+it)/2);
00072     s2m   = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));    // scaling factor for upper range
00073 
00074     overfl  = rbig*s2m;                                             // overflow boundary for abig
00075     eps     = RealScalar(pow(double(ibeta), 1-it));
00076     relerr  = sqrt(eps);                                            // tolerance for neglecting asml
00077     initialized = true;
00078   }
00079   Index n = vec.size();
00080   RealScalar ab2 = b2 / RealScalar(n);
00081   RealScalar asml = RealScalar(0);
00082   RealScalar amed = RealScalar(0);
00083   RealScalar abig = RealScalar(0);
00084   for(typename Derived::InnerIterator it(vec, 0); it; ++it)
00085   {
00086     RealScalar ax = abs(it.value());
00087     if(ax > ab2)     abig += numext::abs2(ax*s2m);
00088     else if(ax < b1) asml += numext::abs2(ax*s1m);
00089     else             amed += numext::abs2(ax);
00090   }
00091   if(abig > RealScalar(0))
00092   {
00093     abig = sqrt(abig);
00094     if(abig > overfl)
00095     {
00096       return rbig;
00097     }
00098     if(amed > RealScalar(0))
00099     {
00100       abig = abig/s2m;
00101       amed = sqrt(amed);
00102     }
00103     else
00104       return abig/s2m;
00105   }
00106   else if(asml > RealScalar(0))
00107   {
00108     if (amed > RealScalar(0))
00109     {
00110       abig = sqrt(amed);
00111       amed = sqrt(asml) / s1m;
00112     }
00113     else
00114       return sqrt(asml)/s1m;
00115   }
00116   else
00117     return sqrt(amed);
00118   asml = (min)(abig, amed);
00119   abig = (max)(abig, amed);
00120   if(asml <= abig*relerr)
00121     return abig;
00122   else
00123     return abig * sqrt(RealScalar(1) + numext::abs2(asml/abig));
00124 }
00125 
00126 } // end namespace internal
00127 
00138 template<typename Derived>
00139 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00140 MatrixBase<Derived>::stableNorm() const
00141 {
00142   using std::min;
00143   using std::sqrt;
00144   const Index blockSize = 4096;
00145   RealScalar scale(0);
00146   RealScalar invScale(1);
00147   RealScalar ssq(0); // sum of square
00148   enum {
00149     Alignment = (int(Flags)&DirectAccessBit) || (int(Flags)&AlignedBit) ? 1 : 0
00150   };
00151   Index n = size();
00152   Index bi = internal::first_aligned(derived());
00153   if (bi>0)
00154     internal::stable_norm_kernel(this->head(bi), ssq, scale, invScale);
00155   for (; bi<n; bi+=blockSize)
00156     internal::stable_norm_kernel(this->segment(bi,(min)(blockSize, n - bi)).template forceAlignedAccessIf<Alignment>(), ssq, scale, invScale);
00157   return scale * sqrt(ssq);
00158 }
00159 
00169 template<typename Derived>
00170 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00171 MatrixBase<Derived>::blueNorm() const
00172 {
00173   return internal::blueNorm_impl(*this);
00174 }
00175 
00181 template<typename Derived>
00182 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
00183 MatrixBase<Derived>::hypotNorm() const
00184 {
00185   return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
00186 }
00187 
00188 } // end namespace Eigen
00189 
00190 #endif // EIGEN_STABLENORM_H


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Thu Aug 27 2015 12:00:57