10 #ifndef EIGEN_STABLENORM_H 11 #define EIGEN_STABLENORM_H 17 template<
typename ExpressionType,
typename Scalar>
18 inline void stable_norm_kernel(
const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
20 Scalar maxCoeff = bl.cwiseAbs().maxCoeff();
24 ssq = ssq * numext::abs2(scale/maxCoeff);
25 Scalar tmp = Scalar(1)/maxCoeff;
26 if(tmp > NumTraits<Scalar>::highest())
28 invScale = NumTraits<Scalar>::highest();
29 scale = Scalar(1)/invScale;
31 else if(maxCoeff>NumTraits<Scalar>::highest())
42 else if(maxCoeff!=maxCoeff)
50 ssq += (bl*invScale).squaredNorm();
53 template<
typename VectorType,
typename RealScalar>
54 void stable_norm_impl_inner_step(
const VectorType &vec, RealScalar& ssq, RealScalar& scale, RealScalar& invScale)
56 typedef typename VectorType::Scalar Scalar;
57 const Index blockSize = 4096;
59 typedef typename internal::nested_eval<VectorType,2>::type VectorTypeCopy;
60 typedef typename internal::remove_all<VectorTypeCopy>::type VectorTypeCopyClean;
61 const VectorTypeCopy copy(vec);
65 || (
int(internal::evaluator<VectorTypeCopyClean>::Alignment)>0)
66 ) && (blockSize*
sizeof(Scalar)*2<EIGEN_STACK_ALLOCATION_LIMIT)
67 && (EIGEN_MAX_STATIC_ALIGN_BYTES>0)
69 typedef typename internal::conditional<CanAlign, Ref<const Matrix<Scalar,Dynamic,1,0,blockSize,1>, internal::evaluator<VectorTypeCopyClean>::Alignment>,
70 typename VectorTypeCopyClean::ConstSegmentReturnType>::type SegmentWrapper;
73 Index bi = internal::first_default_aligned(copy);
75 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, invScale);
80 template<
typename VectorType>
81 typename VectorType::RealScalar
82 stable_norm_impl(
const VectorType &vec,
typename enable_if<VectorType::IsVectorAtCompileTime>::type* = 0 )
90 return abs(vec.coeff(0));
92 typedef typename VectorType::RealScalar RealScalar;
94 RealScalar invScale(1);
97 stable_norm_impl_inner_step(vec, ssq, scale, invScale);
99 return scale *
sqrt(ssq);
102 template<
typename MatrixType>
103 typename MatrixType::RealScalar
104 stable_norm_impl(
const MatrixType &mat,
typename enable_if<!MatrixType::IsVectorAtCompileTime>::type* = 0 )
108 typedef typename MatrixType::RealScalar RealScalar;
110 RealScalar invScale(1);
113 for(
Index j=0; j<mat.outerSize(); ++j)
114 stable_norm_impl_inner_step(mat.innerVector(j), ssq, scale, invScale);
115 return scale *
sqrt(ssq);
118 template<
typename Derived>
119 inline typename NumTraits<typename traits<Derived>::Scalar>::Real
120 blueNorm_impl(
const EigenBase<Derived>& _vec)
122 typedef typename Derived::RealScalar RealScalar;
135 static const int ibeta = std::numeric_limits<RealScalar>::radix;
136 static const int it = NumTraits<RealScalar>::digits();
137 static const int iemin = NumTraits<RealScalar>::min_exponent();
138 static const int iemax = NumTraits<RealScalar>::max_exponent();
139 static const RealScalar rbig = NumTraits<RealScalar>::highest();
140 static const RealScalar b1 = RealScalar(pow(RealScalar(ibeta),RealScalar(-((1-iemin)/2))));
141 static const RealScalar b2 = RealScalar(pow(RealScalar(ibeta),RealScalar((iemax + 1 - it)/2)));
142 static const RealScalar s1m = RealScalar(pow(RealScalar(ibeta),RealScalar((2-iemin)/2)));
143 static const RealScalar s2m = RealScalar(pow(RealScalar(ibeta),RealScalar(- ((iemax+it)/2))));
144 static const RealScalar eps = RealScalar(pow(
double(ibeta), 1-it));
145 static const RealScalar relerr =
sqrt(eps);
147 const Derived& vec(_vec.derived());
148 Index n = vec.size();
149 RealScalar ab2 = b2 / RealScalar(n);
150 RealScalar asml = RealScalar(0);
151 RealScalar amed = RealScalar(0);
152 RealScalar abig = RealScalar(0);
154 for(
Index j=0; j<vec.outerSize(); ++j)
156 for(
typename Derived::InnerIterator iter(vec, j); iter; ++iter)
158 RealScalar ax =
abs(iter.value());
159 if(ax > ab2) abig += numext::abs2(ax*s2m);
160 else if(ax < b1) asml += numext::abs2(ax*s1m);
161 else amed += numext::abs2(ax);
166 if(abig > RealScalar(0))
171 if(amed > RealScalar(0))
179 else if(asml > RealScalar(0))
181 if (amed > RealScalar(0))
184 amed =
sqrt(asml) / s1m;
187 return sqrt(asml)/s1m;
191 asml = numext::mini(abig, amed);
192 abig = numext::maxi(abig, amed);
193 if(asml <= abig*relerr)
196 return abig *
sqrt(RealScalar(1) + numext::abs2(asml/abig));
211 template<
typename Derived>
212 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
215 return internal::stable_norm_impl(derived());
227 template<
typename Derived>
231 return internal::blueNorm_impl(*
this);
239 template<
typename Derived>
244 return numext::abs(coeff(0,0));
246 return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
251 #endif // EIGEN_STABLENORM_H RealScalar blueNorm() const
Definition: StableNorm.h:229
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const unsigned int DirectAccessBit
Definition: Constants.h:155
Namespace containing all symbols from the Eigen library.
Definition: Core:141
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
Definition: NumTraits.h:232
RealScalar stableNorm() const
Definition: StableNorm.h:213
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:74
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
Definition: Eigen_Colamd.h:50
RealScalar hypotNorm() const
Definition: StableNorm.h:241