10 #ifndef EIGEN_STABLENORM_H
11 #define EIGEN_STABLENORM_H
13 #include "./InternalHeaderCheck.h"
19 template<
typename ExpressionType,
typename Scalar>
20 inline void stable_norm_kernel(
const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
22 Scalar maxCoeff = bl.cwiseAbs().maxCoeff();
26 ssq = ssq * numext::abs2(scale/maxCoeff);
27 Scalar tmp = Scalar(1)/maxCoeff;
28 if(tmp > NumTraits<Scalar>::highest())
30 invScale = NumTraits<Scalar>::highest();
31 scale = Scalar(1)/invScale;
33 else if(maxCoeff>NumTraits<Scalar>::highest())
44 else if(maxCoeff!=maxCoeff)
52 ssq += (bl*invScale).squaredNorm();
55 template<
typename VectorType,
typename RealScalar>
56 void stable_norm_impl_inner_step(
const VectorType &vec, RealScalar& ssq, RealScalar& scale, RealScalar& invScale)
58 typedef typename VectorType::Scalar Scalar;
59 const Index blockSize = 4096;
61 typedef typename internal::nested_eval<VectorType,2>::type VectorTypeCopy;
62 typedef internal::remove_all_t<VectorTypeCopy> VectorTypeCopyClean;
63 const VectorTypeCopy copy(vec);
67 || (
int(internal::evaluator<VectorTypeCopyClean>::Alignment)>0)
68 ) && (blockSize*
sizeof(Scalar)*2<EIGEN_STACK_ALLOCATION_LIMIT)
69 && (EIGEN_MAX_STATIC_ALIGN_BYTES>0)
71 typedef std::conditional_t<CanAlign, Ref<const Matrix<Scalar,Dynamic,1,0,blockSize,1>, internal::evaluator<VectorTypeCopyClean>::Alignment>,
72 typename VectorTypeCopyClean::ConstSegmentReturnType> SegmentWrapper;
75 Index bi = internal::first_default_aligned(copy);
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);
82 template<
typename VectorType>
83 typename VectorType::RealScalar
84 stable_norm_impl(
const VectorType &vec, std::enable_if_t<VectorType::IsVectorAtCompileTime>* = 0 )
92 return abs(vec.coeff(0));
94 typedef typename VectorType::RealScalar RealScalar;
96 RealScalar invScale(1);
99 stable_norm_impl_inner_step(vec, ssq, scale, invScale);
101 return scale *
sqrt(ssq);
104 template<
typename MatrixType>
105 typename MatrixType::RealScalar
106 stable_norm_impl(
const MatrixType &mat, std::enable_if_t<!MatrixType::IsVectorAtCompileTime>* = 0 )
110 typedef typename MatrixType::RealScalar RealScalar;
112 RealScalar invScale(1);
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);
120 template<
typename Derived>
121 inline typename NumTraits<typename traits<Derived>::Scalar>::Real
122 blueNorm_impl(
const EigenBase<Derived>& _vec)
124 typedef typename Derived::RealScalar RealScalar;
137 static const int ibeta = std::numeric_limits<RealScalar>::radix;
138 static const int it = NumTraits<RealScalar>::digits();
139 static const int iemin = NumTraits<RealScalar>::min_exponent();
140 static const int iemax = NumTraits<RealScalar>::max_exponent();
141 static const RealScalar rbig = NumTraits<RealScalar>::highest();
142 static const RealScalar b1 = RealScalar(pow(RealScalar(ibeta),RealScalar(-((1-iemin)/2))));
143 static const RealScalar b2 = RealScalar(pow(RealScalar(ibeta),RealScalar((iemax + 1 - it)/2)));
144 static const RealScalar s1m = RealScalar(pow(RealScalar(ibeta),RealScalar((2-iemin)/2)));
145 static const RealScalar s2m = RealScalar(pow(RealScalar(ibeta),RealScalar(- ((iemax+it)/2))));
146 static const RealScalar eps = RealScalar(pow(
double(ibeta), 1-it));
147 static const RealScalar relerr =
sqrt(eps);
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);
156 for(
Index j=0; j<vec.outerSize(); ++j)
158 for(
typename Derived::InnerIterator iter(vec, j); iter; ++iter)
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);
168 if(abig > RealScalar(0))
173 if(amed > RealScalar(0))
181 else if(asml > RealScalar(0))
183 if (amed > RealScalar(0))
186 amed =
sqrt(asml) / s1m;
189 return sqrt(asml)/s1m;
193 asml = numext::mini(abig, amed);
194 abig = numext::maxi(abig, amed);
195 if(asml <= abig*relerr)
198 return abig *
sqrt(RealScalar(1) + numext::abs2(asml/abig));
213 template<
typename Derived>
214 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
217 return internal::stable_norm_impl(derived());
229 template<
typename Derived>
233 return internal::blueNorm_impl(*
this);
241 template<
typename Derived>
246 return numext::abs(coeff(0,0));
248 return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
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
Namespace containing all symbols from the Eigen library.
Definition: Core:139
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:59
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)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:231