StableNorm.h 7.51 KB
Newer Older
LM's avatar
LM committed
1 2 3 4 5
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
//
Don Gagne's avatar
Don Gagne committed
6 7 8
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
LM's avatar
LM committed
9 10 11 12

#ifndef EIGEN_STABLENORM_H
#define EIGEN_STABLENORM_H

Don Gagne's avatar
Don Gagne committed
13 14
namespace Eigen { 

LM's avatar
LM committed
15
namespace internal {
Don Gagne's avatar
Don Gagne committed
16

LM's avatar
LM committed
17 18 19
template<typename ExpressionType, typename Scalar>
inline void stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
{
20 21
  Scalar maxCoeff = bl.cwiseAbs().maxCoeff();
  
22
  if(maxCoeff>scale)
LM's avatar
LM committed
23
  {
24 25 26 27 28 29 30
    ssq = ssq * numext::abs2(scale/maxCoeff);
    Scalar tmp = Scalar(1)/maxCoeff;
    if(tmp > NumTraits<Scalar>::highest())
    {
      invScale = NumTraits<Scalar>::highest();
      scale = Scalar(1)/invScale;
    }
31 32 33 34 35
    else if(maxCoeff>NumTraits<Scalar>::highest()) // we got a INF
    {
      invScale = Scalar(1);
      scale = maxCoeff;
    }
36 37 38 39 40
    else
    {
      scale = maxCoeff;
      invScale = tmp;
    }
LM's avatar
LM committed
41
  }
42 43 44 45
  else if(maxCoeff!=maxCoeff) // we got a NaN
  {
    scale = maxCoeff;
  }
46 47
  
  // TODO if the maxCoeff is much much smaller than the current scale,
LM's avatar
LM committed
48
  // then we can neglect this sub vector
49 50
  if(scale>Scalar(0)) // if scale==0, then bl is 0 
    ssq += (bl*invScale).squaredNorm();
LM's avatar
LM committed
51 52 53
}

template<typename Derived>
Don Gagne's avatar
Don Gagne committed
54 55
inline typename NumTraits<typename traits<Derived>::Scalar>::Real
blueNorm_impl(const EigenBase<Derived>& _vec)
LM's avatar
LM committed
56
{
Don Gagne's avatar
Don Gagne committed
57
  typedef typename Derived::RealScalar RealScalar;  
LM's avatar
LM committed
58
  using std::pow;
Don Gagne's avatar
Don Gagne committed
59 60 61 62
  using std::sqrt;
  using std::abs;
  const Derived& vec(_vec.derived());
  static bool initialized = false;
63
  static RealScalar b1, b2, s1m, s2m, rbig, relerr;
Don Gagne's avatar
Don Gagne committed
64
  if(!initialized)
LM's avatar
LM committed
65
  {
Don Gagne's avatar
Don Gagne committed
66 67
    int ibeta, it, iemin, iemax, iexp;
    RealScalar eps;
LM's avatar
LM committed
68
    // This program calculates the machine-dependent constants
Don Gagne's avatar
Don Gagne committed
69
    // bl, b2, slm, s2m, relerr overfl
LM's avatar
LM committed
70 71 72 73 74 75
    // from the "basic" machine-dependent numbers
    // nbig, ibeta, it, iemin, iemax, rbig.
    // The following define the basic machine-dependent constants.
    // For portability, the PORT subprograms "ilmaeh" and "rlmach"
    // are used. For any specific computer, each of the assignment
    // statements can be replaced
Don Gagne's avatar
Don Gagne committed
76 77 78 79 80
    ibeta = std::numeric_limits<RealScalar>::radix;                 // base for floating-point numbers
    it    = std::numeric_limits<RealScalar>::digits;                // number of base-beta digits in mantissa
    iemin = std::numeric_limits<RealScalar>::min_exponent;          // minimum exponent
    iemax = std::numeric_limits<RealScalar>::max_exponent;          // maximum exponent
    rbig  = (std::numeric_limits<RealScalar>::max)();               // largest floating-point number
LM's avatar
LM committed
81 82

    iexp  = -((1-iemin)/2);
Don Gagne's avatar
Don Gagne committed
83
    b1    = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));    // lower boundary of midrange
LM's avatar
LM committed
84
    iexp  = (iemax + 1 - it)/2;
Don Gagne's avatar
Don Gagne committed
85
    b2    = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));    // upper boundary of midrange
LM's avatar
LM committed
86 87

    iexp  = (2-iemin)/2;
Don Gagne's avatar
Don Gagne committed
88
    s1m   = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));    // scaling factor for lower range
LM's avatar
LM committed
89
    iexp  = - ((iemax+it)/2);
Don Gagne's avatar
Don Gagne committed
90
    s2m   = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp)));    // scaling factor for upper range
LM's avatar
LM committed
91 92

    eps     = RealScalar(pow(double(ibeta), 1-it));
Don Gagne's avatar
Don Gagne committed
93 94
    relerr  = sqrt(eps);                                            // tolerance for neglecting asml
    initialized = true;
LM's avatar
LM committed
95
  }
Don Gagne's avatar
Don Gagne committed
96
  Index n = vec.size();
LM's avatar
LM committed
97 98 99 100
  RealScalar ab2 = b2 / RealScalar(n);
  RealScalar asml = RealScalar(0);
  RealScalar amed = RealScalar(0);
  RealScalar abig = RealScalar(0);
Don Gagne's avatar
Don Gagne committed
101
  for(typename Derived::InnerIterator it(vec, 0); it; ++it)
LM's avatar
LM committed
102
  {
Don Gagne's avatar
Don Gagne committed
103 104 105 106
    RealScalar ax = abs(it.value());
    if(ax > ab2)     abig += numext::abs2(ax*s2m);
    else if(ax < b1) asml += numext::abs2(ax*s1m);
    else             amed += numext::abs2(ax);
LM's avatar
LM committed
107
  }
108 109
  if(amed!=amed)
    return amed;  // we got a NaN
LM's avatar
LM committed
110 111
  if(abig > RealScalar(0))
  {
Don Gagne's avatar
Don Gagne committed
112
    abig = sqrt(abig);
113 114
    if(abig > rbig) // overflow, or *this contains INF values
      return abig;  // return INF
LM's avatar
LM committed
115 116 117
    if(amed > RealScalar(0))
    {
      abig = abig/s2m;
Don Gagne's avatar
Don Gagne committed
118
      amed = sqrt(amed);
LM's avatar
LM committed
119 120 121 122 123 124 125 126
    }
    else
      return abig/s2m;
  }
  else if(asml > RealScalar(0))
  {
    if (amed > RealScalar(0))
    {
Don Gagne's avatar
Don Gagne committed
127 128
      abig = sqrt(amed);
      amed = sqrt(asml) / s1m;
LM's avatar
LM committed
129 130
    }
    else
Don Gagne's avatar
Don Gagne committed
131
      return sqrt(asml)/s1m;
LM's avatar
LM committed
132 133
  }
  else
Don Gagne's avatar
Don Gagne committed
134
    return sqrt(amed);
135 136
  asml = numext::mini(abig, amed);
  abig = numext::maxi(abig, amed);
LM's avatar
LM committed
137 138 139
  if(asml <= abig*relerr)
    return abig;
  else
Don Gagne's avatar
Don Gagne committed
140 141 142 143 144 145 146 147 148 149 150 151 152 153 154 155 156 157 158 159
    return abig * sqrt(RealScalar(1) + numext::abs2(asml/abig));
}

} // end namespace internal

/** \returns the \em l2 norm of \c *this avoiding underflow and overflow.
  * This version use a blockwise two passes algorithm:
  *  1 - find the absolute largest coefficient \c s
  *  2 - compute \f$ s \Vert \frac{*this}{s} \Vert \f$ in a standard way
  *
  * For architecture/scalar types supporting vectorization, this version
  * is faster than blueNorm(). Otherwise the blueNorm() is much faster.
  *
  * \sa norm(), blueNorm(), hypotNorm()
  */
template<typename Derived>
inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
MatrixBase<Derived>::stableNorm() const
{
  using std::sqrt;
160
  using std::abs;
Don Gagne's avatar
Don Gagne committed
161 162 163 164
  const Index blockSize = 4096;
  RealScalar scale(0);
  RealScalar invScale(1);
  RealScalar ssq(0); // sum of square
165 166 167 168 169
  
  typedef typename internal::nested_eval<Derived,2>::type DerivedCopy;
  typedef typename internal::remove_all<DerivedCopy>::type DerivedCopyClean;
  const DerivedCopy copy(derived());
  
Don Gagne's avatar
Don Gagne committed
170
  enum {
171 172 173 174
    CanAlign = (   (int(DerivedCopyClean::Flags)&DirectAccessBit)
                || (int(internal::evaluator<DerivedCopyClean>::Alignment)>0) // FIXME Alignment)>0 might not be enough
               ) && (blockSize*sizeof(Scalar)*2<EIGEN_STACK_ALLOCATION_LIMIT)
                 && (EIGEN_MAX_STATIC_ALIGN_BYTES>0) // if we cannot allocate on the stack, then let's not bother about this optimization
Don Gagne's avatar
Don Gagne committed
175
  };
176 177
  typedef typename internal::conditional<CanAlign, Ref<const Matrix<Scalar,Dynamic,1,0,blockSize,1>, internal::evaluator<DerivedCopyClean>::Alignment>,
                                                   typename DerivedCopyClean::ConstSegmentReturnType>::type SegmentWrapper;
Don Gagne's avatar
Don Gagne committed
178
  Index n = size();
179 180 181 182 183
  
  if(n==1)
    return abs(this->coeff(0));
  
  Index bi = internal::first_default_aligned(copy);
Don Gagne's avatar
Don Gagne committed
184
  if (bi>0)
185
    internal::stable_norm_kernel(copy.head(bi), ssq, scale, invScale);
Don Gagne's avatar
Don Gagne committed
186
  for (; bi<n; bi+=blockSize)
187
    internal::stable_norm_kernel(SegmentWrapper(copy.segment(bi,numext::mini(blockSize, n - bi))), ssq, scale, invScale);
Don Gagne's avatar
Don Gagne committed
188 189 190 191 192 193 194 195 196 197 198 199 200 201 202 203 204
  return scale * sqrt(ssq);
}

/** \returns the \em l2 norm of \c *this using the Blue's algorithm.
  * A Portable Fortran Program to Find the Euclidean Norm of a Vector,
  * ACM TOMS, Vol 4, Issue 1, 1978.
  *
  * For architecture/scalar types without vectorization, this version
  * is much faster than stableNorm(). Otherwise the stableNorm() is faster.
  *
  * \sa norm(), stableNorm(), hypotNorm()
  */
template<typename Derived>
inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
MatrixBase<Derived>::blueNorm() const
{
  return internal::blueNorm_impl(*this);
LM's avatar
LM committed
205 206 207 208 209 210 211 212 213 214 215 216 217 218
}

/** \returns the \em l2 norm of \c *this avoiding undeflow and overflow.
  * This version use a concatenation of hypot() calls, and it is very slow.
  *
  * \sa norm(), stableNorm()
  */
template<typename Derived>
inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
MatrixBase<Derived>::hypotNorm() const
{
  return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
}

Don Gagne's avatar
Don Gagne committed
219 220
} // end namespace Eigen

LM's avatar
LM committed
221
#endif // EIGEN_STABLENORM_H