StableNorm.h
Go to the documentation of this file.
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 #include "./InternalHeaderCheck.h"
14 
15 namespace Eigen {
16 
17 namespace internal {
18 
19 template<typename ExpressionType, typename Scalar>
20 inline void stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
21 {
22  Scalar maxCoeff = bl.cwiseAbs().maxCoeff();
23 
24  if(maxCoeff>scale)
25  {
26  ssq = ssq * numext::abs2(scale/maxCoeff);
27  Scalar tmp = Scalar(1)/maxCoeff;
28  if(tmp > NumTraits<Scalar>::highest())
29  {
30  invScale = NumTraits<Scalar>::highest();
31  scale = Scalar(1)/invScale;
32  }
33  else if(maxCoeff>NumTraits<Scalar>::highest()) // we got a INF
34  {
35  invScale = Scalar(1);
36  scale = maxCoeff;
37  }
38  else
39  {
40  scale = maxCoeff;
41  invScale = tmp;
42  }
43  }
44  else if(maxCoeff!=maxCoeff) // we got a NaN
45  {
46  scale = maxCoeff;
47  }
48 
49  // TODO if the maxCoeff is much much smaller than the current scale,
50  // then we can neglect this sub vector
51  if(scale>Scalar(0)) // if scale==0, then bl is 0
52  ssq += (bl*invScale).squaredNorm();
53 }
54 
55 template<typename VectorType, typename RealScalar>
56 void stable_norm_impl_inner_step(const VectorType &vec, RealScalar& ssq, RealScalar& scale, RealScalar& invScale)
57 {
58  typedef typename VectorType::Scalar Scalar;
59  const Index blockSize = 4096;
60 
61  typedef typename internal::nested_eval<VectorType,2>::type VectorTypeCopy;
62  typedef internal::remove_all_t<VectorTypeCopy> VectorTypeCopyClean;
63  const VectorTypeCopy copy(vec);
64 
65  enum {
66  CanAlign = ( (int(VectorTypeCopyClean::Flags)&DirectAccessBit)
67  || (int(internal::evaluator<VectorTypeCopyClean>::Alignment)>0) // FIXME Alignment)>0 might not be enough
68  ) && (blockSize*sizeof(Scalar)*2<EIGEN_STACK_ALLOCATION_LIMIT)
69  && (EIGEN_MAX_STATIC_ALIGN_BYTES>0) // if we cannot allocate on the stack, then let's not bother about this optimization
70  };
71  typedef std::conditional_t<CanAlign, Ref<const Matrix<Scalar,Dynamic,1,0,blockSize,1>, internal::evaluator<VectorTypeCopyClean>::Alignment>,
72  typename VectorTypeCopyClean::ConstSegmentReturnType> SegmentWrapper;
73  Index n = vec.size();
74 
76  if (bi>0)
77  internal::stable_norm_kernel(copy.head(bi), ssq, scale, invScale);
78  for (; bi<n; bi+=blockSize)
79  internal::stable_norm_kernel(SegmentWrapper(copy.segment(bi,numext::mini(blockSize, n - bi))), ssq, scale, invScale);
80 }
81 
82 template<typename VectorType>
83 typename VectorType::RealScalar
84 stable_norm_impl(const VectorType &vec, std::enable_if_t<VectorType::IsVectorAtCompileTime>* = 0 )
85 {
86  using std::sqrt;
87  using std::abs;
88 
89  Index n = vec.size();
90 
91  if(n==1)
92  return abs(vec.coeff(0));
93 
94  typedef typename VectorType::RealScalar RealScalar;
95  RealScalar scale(0);
96  RealScalar invScale(1);
97  RealScalar ssq(0); // sum of squares
98 
99  stable_norm_impl_inner_step(vec, ssq, scale, invScale);
100 
101  return scale * sqrt(ssq);
102 }
103 
104 template<typename MatrixType>
105 typename MatrixType::RealScalar
106 stable_norm_impl(const MatrixType &mat, std::enable_if_t<!MatrixType::IsVectorAtCompileTime>* = 0 )
107 {
108  using std::sqrt;
109 
110  typedef typename MatrixType::RealScalar RealScalar;
111  RealScalar scale(0);
112  RealScalar invScale(1);
113  RealScalar ssq(0); // sum of squares
114 
115  for(Index j=0; j<mat.outerSize(); ++j)
116  stable_norm_impl_inner_step(mat.innerVector(j), ssq, scale, invScale);
117  return scale * sqrt(ssq);
118 }
119 
120 template<typename Derived>
121 inline typename NumTraits<typename traits<Derived>::Scalar>::Real
123 {
124  typedef typename Derived::RealScalar RealScalar;
125  using std::pow;
126  using std::sqrt;
127  using std::abs;
128 
129  // This program calculates the machine-dependent constants
130  // bl, b2, slm, s2m, relerr overfl
131  // from the "basic" machine-dependent numbers
132  // nbig, ibeta, it, iemin, iemax, rbig.
133  // The following define the basic machine-dependent constants.
134  // For portability, the PORT subprograms "ilmaeh" and "rlmach"
135  // are used. For any specific computer, each of the assignment
136  // statements can be replaced
137  static const int ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers
138  static const int it = NumTraits<RealScalar>::digits(); // number of base-beta digits in mantissa
139  static const int iemin = NumTraits<RealScalar>::min_exponent(); // minimum exponent
140  static const int iemax = NumTraits<RealScalar>::max_exponent(); // maximum exponent
141  static const RealScalar rbig = NumTraits<RealScalar>::highest(); // largest floating-point number
142  static const RealScalar b1 = RealScalar(pow(RealScalar(ibeta),RealScalar(-((1-iemin)/2)))); // lower boundary of midrange
143  static const RealScalar b2 = RealScalar(pow(RealScalar(ibeta),RealScalar((iemax + 1 - it)/2))); // upper boundary of midrange
144  static const RealScalar s1m = RealScalar(pow(RealScalar(ibeta),RealScalar((2-iemin)/2))); // scaling factor for lower range
145  static const RealScalar s2m = RealScalar(pow(RealScalar(ibeta),RealScalar(- ((iemax+it)/2)))); // scaling factor for upper range
146  static const RealScalar eps = RealScalar(pow(double(ibeta), 1-it));
147  static const RealScalar relerr = sqrt(eps); // tolerance for neglecting asml
148 
149  const Derived& vec(_vec.derived());
150  Index n = vec.size();
151  RealScalar ab2 = b2 / RealScalar(n);
152  RealScalar asml = RealScalar(0);
153  RealScalar amed = RealScalar(0);
154  RealScalar abig = RealScalar(0);
155 
156  for(Index j=0; j<vec.outerSize(); ++j)
157  {
158  for(typename Derived::InnerIterator iter(vec, j); iter; ++iter)
159  {
160  RealScalar ax = abs(iter.value());
161  if(ax > ab2) abig += numext::abs2(ax*s2m);
162  else if(ax < b1) asml += numext::abs2(ax*s1m);
163  else amed += numext::abs2(ax);
164  }
165  }
166  if(amed!=amed)
167  return amed; // we got a NaN
168  if(abig > RealScalar(0))
169  {
170  abig = sqrt(abig);
171  if(abig > rbig) // overflow, or *this contains INF values
172  return abig; // return INF
173  if(amed > RealScalar(0))
174  {
175  abig = abig/s2m;
176  amed = sqrt(amed);
177  }
178  else
179  return abig/s2m;
180  }
181  else if(asml > RealScalar(0))
182  {
183  if (amed > RealScalar(0))
184  {
185  abig = sqrt(amed);
186  amed = sqrt(asml) / s1m;
187  }
188  else
189  return sqrt(asml)/s1m;
190  }
191  else
192  return sqrt(amed);
193  asml = numext::mini(abig, amed);
194  abig = numext::maxi(abig, amed);
195  if(asml <= abig*relerr)
196  return abig;
197  else
198  return abig * sqrt(RealScalar(1) + numext::abs2(asml/abig));
199 }
200 
201 } // end namespace internal
202 
213 template<typename Derived>
216 {
217  return internal::stable_norm_impl(derived());
218 }
219 
229 template<typename Derived>
232 {
233  return internal::blueNorm_impl(*this);
234 }
235 
241 template<typename Derived>
244 {
245  if(size()==1)
246  return numext::abs(coeff(0,0));
247  else
248  return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
249 }
250 
251 } // end namespace Eigen
252 
253 #endif // EIGEN_STABLENORM_H
const AbsReturnType abs() const
const SqrtReturnType sqrt() const
int n
#define EIGEN_MAX_STATIC_ALIGN_BYTES
#define EIGEN_STACK_ALLOCATION_LIMIT
Definition: Macros.h:55
const CwiseAbsReturnType cwiseAbs() const
Matrix< float, 1, Dynamic > MatrixType
NumTraits< Scalar >::Real RealScalar
Definition: DenseBase.h:68
RealScalar hypotNorm() const
Definition: StableNorm.h:243
RealScalar blueNorm() const
Definition: StableNorm.h:231
RealScalar stableNorm() const
Definition: StableNorm.h:215
const unsigned int DirectAccessBit
Definition: Constants.h:157
bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
Definition: BFloat16.h:626
VectorType::RealScalar stable_norm_impl(const VectorType &vec, std::enable_if_t< VectorType::IsVectorAtCompileTime > *=0)
Definition: StableNorm.h:84
void stable_norm_impl_inner_step(const VectorType &vec, RealScalar &ssq, RealScalar &scale, RealScalar &invScale)
Definition: StableNorm.h:56
static Index first_default_aligned(const DenseBase< Derived > &m)
typename remove_all< T >::type remove_all_t
Definition: Meta.h:119
NumTraits< typename traits< Derived >::Scalar >::Real blueNorm_impl(const EigenBase< Derived > &_vec)
Definition: StableNorm.h:122
void stable_norm_kernel(const ExpressionType &bl, Scalar &ssq, Scalar &scale, Scalar &invScale)
Definition: StableNorm.h:20
bool abs2(bool x)
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
: InteropHeaders
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:82
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
Derived & derived()
Definition: EigenBase.h:48
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231
std::ptrdiff_t j