StableNorm.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #ifndef EIGEN_STABLENORM_H
11 #define EIGEN_STABLENORM_H
12 
13 namespace Eigen {
14 
15 namespace internal {
16 
17 template<typename ExpressionType, typename Scalar>
18 inline void stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
19 {
20  Scalar max = bl.cwiseAbs().maxCoeff();
21  if (max>scale)
22  {
23  ssq = ssq * abs2(scale/max);
24  scale = max;
25  invScale = Scalar(1)/scale;
26  }
27  // TODO if the max is much much smaller than the current scale,
28  // then we can neglect this sub vector
29  ssq += (bl*invScale).squaredNorm();
30 }
31 }
32 
43 template<typename Derived>
44 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
46 {
47  using std::min;
48  const Index blockSize = 4096;
49  RealScalar scale(0);
50  RealScalar invScale(1);
51  RealScalar ssq(0); // sum of square
52  enum {
53  Alignment = (int(Flags)&DirectAccessBit) || (int(Flags)&AlignedBit) ? 1 : 0
54  };
55  Index n = size();
56  Index bi = internal::first_aligned(derived());
57  if (bi>0)
58  internal::stable_norm_kernel(this->head(bi), ssq, scale, invScale);
59  for (; bi<n; bi+=blockSize)
60  internal::stable_norm_kernel(this->segment(bi,(min)(blockSize, n - bi)).template forceAlignedAccessIf<Alignment>(), ssq, scale, invScale);
61  return scale * internal::sqrt(ssq);
62 }
63 
73 template<typename Derived>
76 {
77  using std::pow;
78  using std::min;
79  using std::max;
80  static bool initialized = false;
81  static RealScalar b1, b2, s1m, s2m, overfl, rbig, relerr;
82  if(!initialized)
83  {
84  int ibeta, it, iemin, iemax, iexp;
85  RealScalar abig, eps;
86  // This program calculates the machine-dependent constants
87  // bl, b2, slm, s2m, relerr overfl
88  // from the "basic" machine-dependent numbers
89  // ibeta, it, iemin, iemax, rbig.
90  // The following define the basic machine-dependent constants.
91  // For portability, the PORT subprograms "ilmaeh" and "rlmach"
92  // are used. For any specific computer, each of the assignment
93  // statements can be replaced
94  ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers
95  it = std::numeric_limits<RealScalar>::digits; // number of base-beta digits in mantissa
96  iemin = std::numeric_limits<RealScalar>::min_exponent; // minimum exponent
97  iemax = std::numeric_limits<RealScalar>::max_exponent; // maximum exponent
98  rbig = (std::numeric_limits<RealScalar>::max)(); // largest floating-point number
99 
100  iexp = -((1-iemin)/2);
101  b1 = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp))); // lower boundary of midrange
102  iexp = (iemax + 1 - it)/2;
103  b2 = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp))); // upper boundary of midrange
104 
105  iexp = (2-iemin)/2;
106  s1m = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp))); // scaling factor for lower range
107  iexp = - ((iemax+it)/2);
108  s2m = RealScalar(pow(RealScalar(ibeta),RealScalar(iexp))); // scaling factor for upper range
109 
110  overfl = rbig*s2m; // overflow boundary for abig
111  eps = RealScalar(pow(double(ibeta), 1-it));
112  relerr = internal::sqrt(eps); // tolerance for neglecting asml
113  abig = RealScalar(1.0/eps - 1.0);
114  initialized = true;
115  }
116  Index n = size();
117  RealScalar ab2 = b2 / RealScalar(n);
118  RealScalar asml = RealScalar(0);
119  RealScalar amed = RealScalar(0);
120  RealScalar abig = RealScalar(0);
121  for(Index j=0; j<n; ++j)
122  {
123  RealScalar ax = internal::abs(coeff(j));
124  if(ax > ab2) abig += internal::abs2(ax*s2m);
125  else if(ax < b1) asml += internal::abs2(ax*s1m);
126  else amed += internal::abs2(ax);
127  }
128  if(abig > RealScalar(0))
129  {
130  abig = internal::sqrt(abig);
131  if(abig > overfl)
132  {
133  return rbig;
134  }
135  if(amed > RealScalar(0))
136  {
137  abig = abig/s2m;
138  amed = internal::sqrt(amed);
139  }
140  else
141  return abig/s2m;
142  }
143  else if(asml > RealScalar(0))
144  {
145  if (amed > RealScalar(0))
146  {
147  abig = internal::sqrt(amed);
148  amed = internal::sqrt(asml) / s1m;
149  }
150  else
151  return internal::sqrt(asml)/s1m;
152  }
153  else
154  return internal::sqrt(amed);
155  asml = (min)(abig, amed);
156  abig = (max)(abig, amed);
157  if(asml <= abig*relerr)
158  return abig;
159  else
160  return abig * internal::sqrt(RealScalar(1) + internal::abs2(asml/abig));
161 }
162 
168 template<typename Derived>
171 {
172  return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
173 }
174 
175 } // end namespace Eigen
176 
177 #endif // EIGEN_STABLENORM_H