Eigen  5.0.1-dev+60122df6
 
Loading...
Searching...
No Matches
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// IWYU pragma: private
14#include "./InternalHeaderCheck.h"
15
16namespace Eigen {
17
18namespace internal {
19
20template <typename ExpressionType, typename Scalar>
21inline void stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale) {
22 Scalar maxCoeff = bl.cwiseAbs().maxCoeff();
23
24 if (maxCoeff > scale) {
25 ssq = ssq * numext::abs2(scale / maxCoeff);
26 Scalar tmp = Scalar(1) / maxCoeff;
27 if (tmp > NumTraits<Scalar>::highest()) {
28 invScale = NumTraits<Scalar>::highest();
29 scale = Scalar(1) / invScale;
30 } else if (maxCoeff > NumTraits<Scalar>::highest()) // we got a INF
31 {
32 invScale = Scalar(1);
33 scale = maxCoeff;
34 } else {
35 scale = maxCoeff;
36 invScale = tmp;
37 }
38 } else if (maxCoeff != maxCoeff) // we got a NaN
39 {
40 scale = maxCoeff;
41 }
42
43 // TODO if the maxCoeff is much much smaller than the current scale,
44 // then we can neglect this sub vector
45 if (scale > Scalar(0)) // if scale==0, then bl is 0
46 ssq += (bl * invScale).squaredNorm();
47}
48
49template <typename VectorType, typename RealScalar>
50void stable_norm_impl_inner_step(const VectorType& vec, RealScalar& ssq, RealScalar& scale, RealScalar& invScale) {
51 const Index blockSize = 4096;
52
53 Index n = vec.size();
54 Index blockEnd = numext::round_down(n, blockSize);
55 for (Index i = 0; i < blockEnd; i += blockSize) {
56 internal::stable_norm_kernel(vec.template segment<blockSize>(i), ssq, scale, invScale);
57 }
58 if (n > blockEnd) {
59 internal::stable_norm_kernel(vec.tail(n - blockEnd), ssq, scale, invScale);
60 }
61}
62
63template <typename VectorType>
64typename VectorType::RealScalar stable_norm_impl(const VectorType& vec,
65 std::enable_if_t<VectorType::IsVectorAtCompileTime>* = 0) {
66 using std::abs;
67 using std::sqrt;
68
69 Index n = vec.size();
70 if (EIGEN_PREDICT_FALSE(n == 1)) return abs(vec.coeff(0));
71
72 typedef typename VectorType::RealScalar RealScalar;
73 RealScalar scale(0);
74 RealScalar invScale(1);
75 RealScalar ssq(0); // sum of squares
76
77 stable_norm_impl_inner_step(vec, ssq, scale, invScale);
78
79 return scale * sqrt(ssq);
80}
81
82template <typename MatrixType>
83typename MatrixType::RealScalar stable_norm_impl(const MatrixType& mat,
84 std::enable_if_t<!MatrixType::IsVectorAtCompileTime>* = 0) {
85 using std::sqrt;
86
87 typedef typename MatrixType::RealScalar RealScalar;
88 RealScalar scale(0);
89 RealScalar invScale(1);
90 RealScalar ssq(0); // sum of squares
91
92 for (Index j = 0; j < mat.outerSize(); ++j) stable_norm_impl_inner_step(mat.innerVector(j), ssq, scale, invScale);
93 return scale * sqrt(ssq);
94}
95
96template <typename Derived>
97inline typename NumTraits<typename traits<Derived>::Scalar>::Real blueNorm_impl(const EigenBase<Derived>& _vec) {
98 typedef typename Derived::RealScalar RealScalar;
99 using std::abs;
100 using std::pow;
101 using std::sqrt;
102
103 // This program calculates the machine-dependent constants
104 // bl, b2, slm, s2m, relerr overfl
105 // from the "basic" machine-dependent numbers
106 // nbig, ibeta, it, iemin, iemax, rbig.
107 // The following define the basic machine-dependent constants.
108 // For portability, the PORT subprograms "ilmaeh" and "rlmach"
109 // are used. For any specific computer, each of the assignment
110 // statements can be replaced
111 static const int ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers
112 static const int it = NumTraits<RealScalar>::digits(); // number of base-beta digits in mantissa
113 static const int iemin = NumTraits<RealScalar>::min_exponent(); // minimum exponent
114 static const int iemax = NumTraits<RealScalar>::max_exponent(); // maximum exponent
115 static const RealScalar rbig = NumTraits<RealScalar>::highest(); // largest floating-point number
116 static const RealScalar b1 =
117 RealScalar(pow(RealScalar(ibeta), RealScalar(-((1 - iemin) / 2)))); // lower boundary of midrange
118 static const RealScalar b2 =
119 RealScalar(pow(RealScalar(ibeta), RealScalar((iemax + 1 - it) / 2))); // upper boundary of midrange
120 static const RealScalar s1m =
121 RealScalar(pow(RealScalar(ibeta), RealScalar((2 - iemin) / 2))); // scaling factor for lower range
122 static const RealScalar s2m =
123 RealScalar(pow(RealScalar(ibeta), RealScalar(-((iemax + it) / 2)))); // scaling factor for upper range
124 static const RealScalar eps = RealScalar(pow(double(ibeta), 1 - it));
125 static const RealScalar relerr = sqrt(eps); // tolerance for neglecting asml
126
127 const Derived& vec(_vec.derived());
128 Index n = vec.size();
129 RealScalar ab2 = b2 / RealScalar(n);
130 RealScalar asml = RealScalar(0);
131 RealScalar amed = RealScalar(0);
132 RealScalar abig = RealScalar(0);
133
134 for (Index j = 0; j < vec.outerSize(); ++j) {
135 for (typename Derived::InnerIterator iter(vec, j); iter; ++iter) {
136 RealScalar ax = abs(iter.value());
137 if (ax > ab2)
138 abig += numext::abs2(ax * s2m);
139 else if (ax < b1)
140 asml += numext::abs2(ax * s1m);
141 else
142 amed += numext::abs2(ax);
143 }
144 }
145 if (amed != amed) return amed; // we got a NaN
146 if (abig > RealScalar(0)) {
147 abig = sqrt(abig);
148 if (abig > rbig) // overflow, or *this contains INF values
149 return abig; // return INF
150 if (amed > RealScalar(0)) {
151 abig = abig / s2m;
152 amed = sqrt(amed);
153 } else
154 return abig / s2m;
155 } else if (asml > RealScalar(0)) {
156 if (amed > RealScalar(0)) {
157 abig = sqrt(amed);
158 amed = sqrt(asml) / s1m;
159 } else
160 return sqrt(asml) / s1m;
161 } else
162 return sqrt(amed);
163 asml = numext::mini(abig, amed);
164 abig = numext::maxi(abig, amed);
165 if (asml <= abig * relerr)
166 return abig;
167 else
168 return abig * sqrt(RealScalar(1) + numext::abs2(asml / abig));
169}
170
171} // end namespace internal
172
183template <typename Derived>
185 return internal::stable_norm_impl(derived());
186}
187
189 * A Portable Fortran Program to Find the Euclidean Norm of a Vector,
190 * ACM TOMS, Vol 4, Issue 1, 1978.
192 * For architecture/scalar types without vectorization, this version
193 * is much faster than stableNorm(). Otherwise the stableNorm() is faster.
194 *
195 * \sa norm(), stableNorm(), hypotNorm()
196 */
197template <typename Derived>
199 return internal::blueNorm_impl(*this);
200}
201
207template <typename Derived>
209 if (size() == 1)
210 return numext::abs(coeff(0, 0));
211 else
212 return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
213}
214
215} // end namespace Eigen
216
217#endif // EIGEN_STABLENORM_H
const GlobalUnaryPowReturnType< Derived, ScalarExponent > pow(const Eigen::ArrayBase< Derived > &x, const ScalarExponent &exponent)
RealScalar blueNorm() const
Definition StableNorm.h:198
RealScalar hypotNorm() const
Definition StableNorm.h:208
RealScalar stableNorm() const
Definition StableNorm.h:184
const CwiseAbsReturnType cwiseAbs() const
Definition MatrixBase.h:34
Namespace containing all symbols from the Eigen library.
Definition B01_Experimental.dox:1
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition Meta.h:82
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition NumTraits.h:232