Back to home page

EIC code displayed by LXR

 
 

    


File indexing completed on 2025-01-18 09:56:18

0001 // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 #ifndef EIGEN_STABLENORM_H
0011 #define EIGEN_STABLENORM_H
0012 
0013 namespace Eigen { 
0014 
0015 namespace internal {
0016 
0017 template<typename ExpressionType, typename Scalar>
0018 inline void stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
0019 {
0020   Scalar maxCoeff = bl.cwiseAbs().maxCoeff();
0021   
0022   if(maxCoeff>scale)
0023   {
0024     ssq = ssq * numext::abs2(scale/maxCoeff);
0025     Scalar tmp = Scalar(1)/maxCoeff;
0026     if(tmp > NumTraits<Scalar>::highest())
0027     {
0028       invScale = NumTraits<Scalar>::highest();
0029       scale = Scalar(1)/invScale;
0030     }
0031     else if(maxCoeff>NumTraits<Scalar>::highest()) // we got a INF
0032     {
0033       invScale = Scalar(1);
0034       scale = maxCoeff;
0035     }
0036     else
0037     {
0038       scale = maxCoeff;
0039       invScale = tmp;
0040     }
0041   }
0042   else if(maxCoeff!=maxCoeff) // we got a NaN
0043   {
0044     scale = maxCoeff;
0045   }
0046   
0047   // TODO if the maxCoeff is much much smaller than the current scale,
0048   // then we can neglect this sub vector
0049   if(scale>Scalar(0)) // if scale==0, then bl is 0 
0050     ssq += (bl*invScale).squaredNorm();
0051 }
0052 
0053 template<typename VectorType, typename RealScalar>
0054 void stable_norm_impl_inner_step(const VectorType &vec, RealScalar& ssq, RealScalar& scale, RealScalar& invScale)
0055 {
0056   typedef typename VectorType::Scalar Scalar;
0057   const Index blockSize = 4096;
0058   
0059   typedef typename internal::nested_eval<VectorType,2>::type VectorTypeCopy;
0060   typedef typename internal::remove_all<VectorTypeCopy>::type VectorTypeCopyClean;
0061   const VectorTypeCopy copy(vec);
0062   
0063   enum {
0064     CanAlign = (   (int(VectorTypeCopyClean::Flags)&DirectAccessBit)
0065                 || (int(internal::evaluator<VectorTypeCopyClean>::Alignment)>0) // FIXME Alignment)>0 might not be enough
0066                ) && (blockSize*sizeof(Scalar)*2<EIGEN_STACK_ALLOCATION_LIMIT)
0067                  && (EIGEN_MAX_STATIC_ALIGN_BYTES>0) // if we cannot allocate on the stack, then let's not bother about this optimization
0068   };
0069   typedef typename internal::conditional<CanAlign, Ref<const Matrix<Scalar,Dynamic,1,0,blockSize,1>, internal::evaluator<VectorTypeCopyClean>::Alignment>,
0070                                                    typename VectorTypeCopyClean::ConstSegmentReturnType>::type SegmentWrapper;
0071   Index n = vec.size();
0072   
0073   Index bi = internal::first_default_aligned(copy);
0074   if (bi>0)
0075     internal::stable_norm_kernel(copy.head(bi), ssq, scale, invScale);
0076   for (; bi<n; bi+=blockSize)
0077     internal::stable_norm_kernel(SegmentWrapper(copy.segment(bi,numext::mini(blockSize, n - bi))), ssq, scale, invScale);
0078 }
0079 
0080 template<typename VectorType>
0081 typename VectorType::RealScalar
0082 stable_norm_impl(const VectorType &vec, typename enable_if<VectorType::IsVectorAtCompileTime>::type* = 0 )
0083 {
0084   using std::sqrt;
0085   using std::abs;
0086 
0087   Index n = vec.size();
0088 
0089   if(n==1)
0090     return abs(vec.coeff(0));
0091 
0092   typedef typename VectorType::RealScalar RealScalar;
0093   RealScalar scale(0);
0094   RealScalar invScale(1);
0095   RealScalar ssq(0); // sum of squares
0096 
0097   stable_norm_impl_inner_step(vec, ssq, scale, invScale);
0098   
0099   return scale * sqrt(ssq);
0100 }
0101 
0102 template<typename MatrixType>
0103 typename MatrixType::RealScalar
0104 stable_norm_impl(const MatrixType &mat, typename enable_if<!MatrixType::IsVectorAtCompileTime>::type* = 0 )
0105 {
0106   using std::sqrt;
0107 
0108   typedef typename MatrixType::RealScalar RealScalar;
0109   RealScalar scale(0);
0110   RealScalar invScale(1);
0111   RealScalar ssq(0); // sum of squares
0112 
0113   for(Index j=0; j<mat.outerSize(); ++j)
0114     stable_norm_impl_inner_step(mat.innerVector(j), ssq, scale, invScale);
0115   return scale * sqrt(ssq);
0116 }
0117 
0118 template<typename Derived>
0119 inline typename NumTraits<typename traits<Derived>::Scalar>::Real
0120 blueNorm_impl(const EigenBase<Derived>& _vec)
0121 {
0122   typedef typename Derived::RealScalar RealScalar;  
0123   using std::pow;
0124   using std::sqrt;
0125   using std::abs;
0126 
0127   // This program calculates the machine-dependent constants
0128   // bl, b2, slm, s2m, relerr overfl
0129   // from the "basic" machine-dependent numbers
0130   // nbig, ibeta, it, iemin, iemax, rbig.
0131   // The following define the basic machine-dependent constants.
0132   // For portability, the PORT subprograms "ilmaeh" and "rlmach"
0133   // are used. For any specific computer, each of the assignment
0134   // statements can be replaced
0135   static const int ibeta = std::numeric_limits<RealScalar>::radix;  // base for floating-point numbers
0136   static const int it    = NumTraits<RealScalar>::digits();  // number of base-beta digits in mantissa
0137   static const int iemin = NumTraits<RealScalar>::min_exponent();  // minimum exponent
0138   static const int iemax = NumTraits<RealScalar>::max_exponent();  // maximum exponent
0139   static const RealScalar rbig   = NumTraits<RealScalar>::highest();  // largest floating-point number
0140   static const RealScalar b1     = RealScalar(pow(RealScalar(ibeta),RealScalar(-((1-iemin)/2))));  // lower boundary of midrange
0141   static const RealScalar b2     = RealScalar(pow(RealScalar(ibeta),RealScalar((iemax + 1 - it)/2)));  // upper boundary of midrange
0142   static const RealScalar s1m    = RealScalar(pow(RealScalar(ibeta),RealScalar((2-iemin)/2)));  // scaling factor for lower range
0143   static const RealScalar s2m    = RealScalar(pow(RealScalar(ibeta),RealScalar(- ((iemax+it)/2))));  // scaling factor for upper range
0144   static const RealScalar eps    = RealScalar(pow(double(ibeta), 1-it));
0145   static const RealScalar relerr = sqrt(eps);  // tolerance for neglecting asml
0146 
0147   const Derived& vec(_vec.derived());
0148   Index n = vec.size();
0149   RealScalar ab2 = b2 / RealScalar(n);
0150   RealScalar asml = RealScalar(0);
0151   RealScalar amed = RealScalar(0);
0152   RealScalar abig = RealScalar(0);
0153 
0154   for(Index j=0; j<vec.outerSize(); ++j)
0155   {
0156     for(typename Derived::InnerIterator iter(vec, j); iter; ++iter)
0157     {
0158       RealScalar ax = abs(iter.value());
0159       if(ax > ab2)     abig += numext::abs2(ax*s2m);
0160       else if(ax < b1) asml += numext::abs2(ax*s1m);
0161       else             amed += numext::abs2(ax);
0162     }
0163   }
0164   if(amed!=amed)
0165     return amed;  // we got a NaN
0166   if(abig > RealScalar(0))
0167   {
0168     abig = sqrt(abig);
0169     if(abig > rbig) // overflow, or *this contains INF values
0170       return abig;  // return INF
0171     if(amed > RealScalar(0))
0172     {
0173       abig = abig/s2m;
0174       amed = sqrt(amed);
0175     }
0176     else
0177       return abig/s2m;
0178   }
0179   else if(asml > RealScalar(0))
0180   {
0181     if (amed > RealScalar(0))
0182     {
0183       abig = sqrt(amed);
0184       amed = sqrt(asml) / s1m;
0185     }
0186     else
0187       return sqrt(asml)/s1m;
0188   }
0189   else
0190     return sqrt(amed);
0191   asml = numext::mini(abig, amed);
0192   abig = numext::maxi(abig, amed);
0193   if(asml <= abig*relerr)
0194     return abig;
0195   else
0196     return abig * sqrt(RealScalar(1) + numext::abs2(asml/abig));
0197 }
0198 
0199 } // end namespace internal
0200 
0201 /** \returns the \em l2 norm of \c *this avoiding underflow and overflow.
0202   * This version use a blockwise two passes algorithm:
0203   *  1 - find the absolute largest coefficient \c s
0204   *  2 - compute \f$ s \Vert \frac{*this}{s} \Vert \f$ in a standard way
0205   *
0206   * For architecture/scalar types supporting vectorization, this version
0207   * is faster than blueNorm(). Otherwise the blueNorm() is much faster.
0208   *
0209   * \sa norm(), blueNorm(), hypotNorm()
0210   */
0211 template<typename Derived>
0212 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
0213 MatrixBase<Derived>::stableNorm() const
0214 {
0215   return internal::stable_norm_impl(derived());
0216 }
0217 
0218 /** \returns the \em l2 norm of \c *this using the Blue's algorithm.
0219   * A Portable Fortran Program to Find the Euclidean Norm of a Vector,
0220   * ACM TOMS, Vol 4, Issue 1, 1978.
0221   *
0222   * For architecture/scalar types without vectorization, this version
0223   * is much faster than stableNorm(). Otherwise the stableNorm() is faster.
0224   *
0225   * \sa norm(), stableNorm(), hypotNorm()
0226   */
0227 template<typename Derived>
0228 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
0229 MatrixBase<Derived>::blueNorm() const
0230 {
0231   return internal::blueNorm_impl(*this);
0232 }
0233 
0234 /** \returns the \em l2 norm of \c *this avoiding undeflow and overflow.
0235   * This version use a concatenation of hypot() calls, and it is very slow.
0236   *
0237   * \sa norm(), stableNorm()
0238   */
0239 template<typename Derived>
0240 inline typename NumTraits<typename internal::traits<Derived>::Scalar>::Real
0241 MatrixBase<Derived>::hypotNorm() const
0242 {
0243   if(size()==1)
0244     return numext::abs(coeff(0,0));
0245   else
0246     return this->cwiseAbs().redux(internal::scalar_hypot_op<RealScalar>());
0247 }
0248 
0249 } // end namespace Eigen
0250 
0251 #endif // EIGEN_STABLENORM_H