10 #ifndef EIGEN_STABLENORM_H
11 #define EIGEN_STABLENORM_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();
27 Scalar tmp = Scalar(1)/maxCoeff;
31 scale = Scalar(1)/invScale;
44 else if(maxCoeff!=maxCoeff)
52 ssq += (bl*invScale).squaredNorm();
55 template<
typename VectorType,
typename RealScalar>
58 typedef typename VectorType::Scalar Scalar;
59 const Index blockSize = 4096;
61 typedef typename internal::nested_eval<VectorType,2>::type VectorTypeCopy;
63 const VectorTypeCopy copy(vec);
67 || (
int(internal::evaluator<VectorTypeCopyClean>::Alignment)>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;
78 for (; bi<
n; bi+=blockSize)
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);
101 return scale *
sqrt(ssq);
104 template<
typename MatrixType>
112 RealScalar invScale(1);
117 return scale *
sqrt(ssq);
120 template<
typename Derived>
124 typedef typename Derived::RealScalar RealScalar;
137 static const int ibeta = std::numeric_limits<RealScalar>::radix;
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());
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());
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;
195 if(asml <= abig*relerr)
213 template<
typename Derived>
229 template<
typename Derived>
241 template<
typename Derived>
248 return this->
cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
const AbsReturnType abs() const
const SqrtReturnType sqrt() const
#define EIGEN_STACK_ALLOCATION_LIMIT
const CwiseAbsReturnType cwiseAbs() const
Matrix< float, 1, Dynamic > MatrixType
NumTraits< Scalar >::Real RealScalar
RealScalar hypotNorm() const
RealScalar blueNorm() const
RealScalar stableNorm() const
const unsigned int DirectAccessBit
bfloat16 pow(const bfloat16 &a, const bfloat16 &b)
VectorType::RealScalar stable_norm_impl(const VectorType &vec, std::enable_if_t< VectorType::IsVectorAtCompileTime > *=0)
void stable_norm_impl_inner_step(const VectorType &vec, RealScalar &ssq, RealScalar &scale, RealScalar &invScale)
static Index first_default_aligned(const DenseBase< Derived > &m)
typename remove_all< T >::type remove_all_t
NumTraits< typename traits< Derived >::Scalar >::Real blueNorm_impl(const EigenBase< Derived > &_vec)
void stable_norm_kernel(const ExpressionType &bl, Scalar &ssq, Scalar &scale, Scalar &invScale)
EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
EIGEN_ALWAYS_INLINE std::enable_if_t< NumTraits< T >::IsSigned||NumTraits< T >::IsComplex, typename NumTraits< T >::Real > abs(const T &x)
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
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.