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