Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
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 <[email protected]>
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 typedef typename VectorType::Scalar Scalar;
52 const Index blockSize = 4096;
53
54 typedef typename internal::nested_eval<VectorType, 2>::type VectorTypeCopy;
55 typedef internal::remove_all_t<VectorTypeCopy> VectorTypeCopyClean;
56 const VectorTypeCopy copy(vec);
57
58 enum {
59 CanAlign =
60 ((int(VectorTypeCopyClean::Flags) & DirectAccessBit) ||
61 (int(internal::evaluator<VectorTypeCopyClean>::Alignment) > 0) // FIXME Alignment)>0 might not be enough
62 ) &&
63 (blockSize * sizeof(Scalar) * 2 < EIGEN_STACK_ALLOCATION_LIMIT) &&
64 (EIGEN_MAX_STATIC_ALIGN_BYTES >
65 0) // if we cannot allocate on the stack, then let's not bother about this optimization
66 };
67 typedef std::conditional_t<
68 CanAlign,
69 Ref<const Matrix<Scalar, Dynamic, 1, 0, blockSize, 1>, internal::evaluator<VectorTypeCopyClean>::Alignment>,
70 typename VectorTypeCopyClean::ConstSegmentReturnType>
71 SegmentWrapper;
72 Index n = vec.size();
73
74 Index bi = internal::first_default_aligned(copy);
75 if (bi > 0) internal::stable_norm_kernel(copy.head(bi), ssq, scale, invScale);
76 for (; bi < n; bi += blockSize)
77 internal::stable_norm_kernel(SegmentWrapper(copy.segment(bi, numext::mini(blockSize, n - bi))), ssq, scale,
78 invScale);
79}
80
81template <typename VectorType>
82typename VectorType::RealScalar stable_norm_impl(const VectorType& vec,
83 std::enable_if_t<VectorType::IsVectorAtCompileTime>* = 0) {
84 using std::abs;
85 using std::sqrt;
86
87 Index n = vec.size();
88
89 if (n == 1) return abs(vec.coeff(0));
90
91 typedef typename VectorType::RealScalar RealScalar;
92 RealScalar scale(0);
93 RealScalar invScale(1);
94 RealScalar ssq(0); // sum of squares
95
96 stable_norm_impl_inner_step(vec, ssq, scale, invScale);
97
98 return scale * sqrt(ssq);
99}
100
101template <typename MatrixType>
102typename MatrixType::RealScalar stable_norm_impl(const MatrixType& mat,
103 std::enable_if_t<!MatrixType::IsVectorAtCompileTime>* = 0) {
104 using std::sqrt;
105
106 typedef typename MatrixType::RealScalar RealScalar;
107 RealScalar scale(0);
108 RealScalar invScale(1);
109 RealScalar ssq(0); // sum of squares
110
111 for (Index j = 0; j < mat.outerSize(); ++j) stable_norm_impl_inner_step(mat.innerVector(j), ssq, scale, invScale);
112 return scale * sqrt(ssq);
113}
114
115template <typename Derived>
116inline typename NumTraits<typename traits<Derived>::Scalar>::Real blueNorm_impl(const EigenBase<Derived>& _vec) {
117 typedef typename Derived::RealScalar RealScalar;
118 using std::abs;
119 using std::pow;
120 using std::sqrt;
121
122 // This program calculates the machine-dependent constants
123 // bl, b2, slm, s2m, relerr overfl
124 // from the "basic" machine-dependent numbers
125 // nbig, ibeta, it, iemin, iemax, rbig.
126 // The following define the basic machine-dependent constants.
127 // For portability, the PORT subprograms "ilmaeh" and "rlmach"
128 // are used. For any specific computer, each of the assignment
129 // statements can be replaced
130 static const int ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers
131 static const int it = NumTraits<RealScalar>::digits(); // number of base-beta digits in mantissa
132 static const int iemin = NumTraits<RealScalar>::min_exponent(); // minimum exponent
133 static const int iemax = NumTraits<RealScalar>::max_exponent(); // maximum exponent
134 static const RealScalar rbig = NumTraits<RealScalar>::highest(); // largest floating-point number
135 static const RealScalar b1 =
136 RealScalar(pow(RealScalar(ibeta), RealScalar(-((1 - iemin) / 2)))); // lower boundary of midrange
137 static const RealScalar b2 =
138 RealScalar(pow(RealScalar(ibeta), RealScalar((iemax + 1 - it) / 2))); // upper boundary of midrange
139 static const RealScalar s1m =
140 RealScalar(pow(RealScalar(ibeta), RealScalar((2 - iemin) / 2))); // scaling factor for lower range
141 static const RealScalar s2m =
142 RealScalar(pow(RealScalar(ibeta), RealScalar(-((iemax + it) / 2)))); // scaling factor for upper range
143 static const RealScalar eps = RealScalar(pow(double(ibeta), 1 - it));
144 static const RealScalar relerr = sqrt(eps); // tolerance for neglecting asml
145
146 const Derived& vec(_vec.derived());
147 Index n = vec.size();
148 RealScalar ab2 = b2 / RealScalar(n);
149 RealScalar asml = RealScalar(0);
150 RealScalar amed = RealScalar(0);
151 RealScalar abig = RealScalar(0);
152
153 for (Index j = 0; j < vec.outerSize(); ++j) {
154 for (typename Derived::InnerIterator iter(vec, j); iter; ++iter) {
155 RealScalar ax = abs(iter.value());
156 if (ax > ab2)
157 abig += numext::abs2(ax * s2m);
158 else if (ax < b1)
159 asml += numext::abs2(ax * s1m);
160 else
161 amed += numext::abs2(ax);
162 }
163 }
164 if (amed != amed) return amed; // we got a NaN
165 if (abig > RealScalar(0)) {
166 abig = sqrt(abig);
167 if (abig > rbig) // overflow, or *this contains INF values
168 return abig; // return INF
169 if (amed > RealScalar(0)) {
170 abig = abig / s2m;
171 amed = sqrt(amed);
172 } else
173 return abig / s2m;
174 } else if (asml > RealScalar(0)) {
175 if (amed > RealScalar(0)) {
176 abig = sqrt(amed);
177 amed = sqrt(asml) / s1m;
178 } else
179 return sqrt(asml) / s1m;
180 } else
181 return sqrt(amed);
182 asml = numext::mini(abig, amed);
183 abig = numext::maxi(abig, amed);
184 if (asml <= abig * relerr)
185 return abig;
186 else
187 return abig * sqrt(RealScalar(1) + numext::abs2(asml / abig));
188}
190} // end namespace internal
202template <typename Derived>
204 return internal::stable_norm_impl(derived());
205}
206
216template <typename Derived>
218 return internal::blueNorm_impl(*this);
219}
220
226template <typename Derived>
228 if (size() == 1)
229 return numext::abs(coeff(0, 0));
230 else
231 return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
232}
233
234} // end namespace Eigen
235
236#endif // EIGEN_STABLENORM_H
Base class for all dense matrices, vectors, and expressions.
Definition MatrixBase.h:52
const unsigned int DirectAccessBit
Definition Constants.h:159
Namespace containing all symbols from the Eigen library.
Definition Core:137
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)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition Meta.h:523