aboutsummaryrefslogtreecommitdiffhomepage
path: root/Eigen/src/Core
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2009-07-28 17:35:07 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2009-07-28 17:35:07 +0200
commit54804eb62642ab1be510e41db9b573c6f6151bf2 (patch)
tree76eda2dedb4a66be0072425ebd110546211f1f71 /Eigen/src/Core
parent264fe82c655a26f3c3ab5057684dbc51cf533056 (diff)
parent562864bcfb363f603f40ce716c49539fcd1565d3 (diff)
synch with main branch
Diffstat (limited to 'Eigen/src/Core')
-rw-r--r--Eigen/src/Core/Block.h4
-rw-r--r--Eigen/src/Core/Dot.h148
-rw-r--r--Eigen/src/Core/Functors.h3
-rw-r--r--Eigen/src/Core/MapBase.h10
-rw-r--r--Eigen/src/Core/Matrix.h10
-rw-r--r--Eigen/src/Core/MatrixBase.h1
-rw-r--r--Eigen/src/Core/NumTraits.h8
-rw-r--r--Eigen/src/Core/Product.h12
-rw-r--r--Eigen/src/Core/StableNorm.h194
9 files changed, 227 insertions, 163 deletions
diff --git a/Eigen/src/Core/Block.h b/Eigen/src/Core/Block.h
index ffe065e8b..cf7730170 100644
--- a/Eigen/src/Core/Block.h
+++ b/Eigen/src/Core/Block.h
@@ -33,8 +33,8 @@
* \param MatrixType the type of the object in which we are taking a block
* \param BlockRows the number of rows of the block we are taking at compile time (optional)
* \param BlockCols the number of columns of the block we are taking at compile time (optional)
- * \param _PacketAccess allows to enforce aligned loads and stores if set to ForceAligned.
- * The default is AsRequested. This parameter is internaly used by Eigen
+ * \param _PacketAccess allows to enforce aligned loads and stores if set to \b ForceAligned.
+ * The default is \b AsRequested. This parameter is internaly used by Eigen
* in expressions such as \code mat.block() += other; \endcode and most of
* the time this is the only way it is used.
* \param _DirectAccessStatus \internal used for partial specialization
diff --git a/Eigen/src/Core/Dot.h b/Eigen/src/Core/Dot.h
index c5f2e8505..9e84d72bb 100644
--- a/Eigen/src/Core/Dot.h
+++ b/Eigen/src/Core/Dot.h
@@ -292,154 +292,6 @@ inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real MatrixBase<
return ei_sqrt(squaredNorm());
}
-/** \returns the \em l2 norm of \c *this using a numerically more stable
- * algorithm.
- *
- * \sa norm(), dot(), squaredNorm(), blueNorm()
- */
-template<typename Derived>
-inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real
-MatrixBase<Derived>::stableNorm() const
-{
- return this->cwise().abs().redux(ei_scalar_hypot_op<RealScalar>());
-}
-
-/** \internal Computes ibeta^iexp by binary expansion of iexp,
- * exact if ibeta is the machine base */
-template<typename T> inline T bexp(int ibeta, int iexp)
-{
- T tbeta = T(ibeta);
- T res = 1.0;
- int n = iexp;
- if (n<0)
- {
- n = - n;
- tbeta = 1.0/tbeta;
- }
- for(;;)
- {
- if ((n % 2)==0)
- res = res * tbeta;
- n = n/2;
- if (n==0) return res;
- tbeta = tbeta*tbeta;
- }
- return res;
-}
-
-/** \returns the \em l2 norm of \c *this using the Blue's algorithm.
- * A Portable Fortran Program to Find the Euclidean Norm of a Vector,
- * ACM TOMS, Vol 4, Issue 1, 1978.
- *
- * \sa norm(), dot(), squaredNorm(), stableNorm()
- */
-template<typename Derived>
-inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real
-MatrixBase<Derived>::blueNorm() const
-{
- static int nmax;
- static Scalar b1, b2, s1m, s2m, overfl, rbig, relerr;
- int n;
- Scalar ax, abig, amed, asml;
-
- if(nmax <= 0)
- {
- int nbig, ibeta, it, iemin, iemax, iexp;
- Scalar abig, eps;
- // This program calculates the machine-dependent constants
- // bl, b2, slm, s2m, relerr overfl, nmax
- // from the "basic" machine-dependent numbers
- // nbig, ibeta, it, iemin, iemax, rbig.
- // The following define the basic machine-dependent constants.
- // For portability, the PORT subprograms "ilmaeh" and "rlmach"
- // are used. For any specific computer, each of the assignment
- // statements can be replaced
- nbig = std::numeric_limits<int>::max(); // largest integer
- ibeta = NumTraits<Scalar>::Base; // base for floating-point numbers
- it = NumTraits<Scalar>::Mantissa; // number of base-beta digits in mantissa
- iemin = std::numeric_limits<Scalar>::min_exponent; // minimum exponent
- iemax = std::numeric_limits<Scalar>::max_exponent; // maximum exponent
- rbig = std::numeric_limits<Scalar>::max(); // largest floating-point number
-
- // Check the basic machine-dependent constants.
- if(iemin > 1 - 2*it || 1+it>iemax || (it==2 && ibeta<5)
- || (it<=4 && ibeta <= 3 ) || it<2)
- {
- ei_assert(false && "the algorithm cannot be guaranteed on this computer");
- }
- iexp = -((1-iemin)/2);
- b1 = bexp<Scalar>(ibeta, iexp); // lower boundary of midrange
- iexp = (iemax + 1 - it)/2;
- b2 = bexp<Scalar>(ibeta,iexp); // upper boundary of midrange
-
- iexp = (2-iemin)/2;
- s1m = bexp<Scalar>(ibeta,iexp); // scaling factor for lower range
- iexp = - ((iemax+it)/2);
- s2m = bexp<Scalar>(ibeta,iexp); // scaling factor for upper range
-
- overfl = rbig*s2m; // overfow boundary for abig
- eps = bexp<Scalar>(ibeta, 1-it);
- relerr = ei_sqrt(eps); // tolerance for neglecting asml
- abig = 1.0/eps - 1.0;
- if (Scalar(nbig)>abig) nmax = abig; // largest safe n
- else nmax = nbig;
- }
- n = size();
- if(n==0)
- return 0;
- asml = Scalar(0);
- amed = Scalar(0);
- abig = Scalar(0);
- for(int j=0; j<n; ++j)
- {
- ax = ei_abs(coeff(j));
- if(ax > b2) abig += ei_abs2(ax*s2m);
- else if(ax < b1) asml += ei_abs2(ax*s1m);
- else amed += ei_abs2(ax);
- }
- if(abig > Scalar(0))
- {
- abig = ei_sqrt(abig);
- if(abig > overfl)
- {
- ei_assert(false && "overflow");
- return rbig;
- }
- if(amed > Scalar(0))
- {
- abig = abig/s2m;
- amed = ei_sqrt(amed);
- }
- else
- {
- return abig/s2m;
- }
-
- }
- else if(asml > Scalar(0))
- {
- if (amed > Scalar(0))
- {
- abig = ei_sqrt(amed);
- amed = ei_sqrt(asml) / s1m;
- }
- else
- {
- return ei_sqrt(asml)/s1m;
- }
- }
- else
- {
- return ei_sqrt(amed);
- }
- asml = std::min(abig, amed);
- abig = std::max(abig, amed);
- if(asml <= abig*relerr)
- return abig;
- else
- return abig * ei_sqrt(Scalar(1) + ei_abs2(asml/abig));
-}
-
/** \returns an expression of the quotient of *this by its own norm.
*
* \only_for_vectors
diff --git a/Eigen/src/Core/Functors.h b/Eigen/src/Core/Functors.h
index 89badb353..a4c9604df 100644
--- a/Eigen/src/Core/Functors.h
+++ b/Eigen/src/Core/Functors.h
@@ -124,9 +124,6 @@ template<typename Scalar> struct ei_scalar_hypot_op EIGEN_EMPTY_STRUCT {
// typedef typename NumTraits<Scalar>::Real result_type;
EIGEN_STRONG_INLINE const Scalar operator() (const Scalar& _x, const Scalar& _y) const
{
-// typedef typename NumTraits<T>::Real RealScalar;
-// RealScalar _x = ei_abs(x);
-// RealScalar _y = ei_abs(y);
Scalar p = std::max(_x, _y);
Scalar q = std::min(_x, _y);
Scalar qp = q/p;
diff --git a/Eigen/src/Core/MapBase.h b/Eigen/src/Core/MapBase.h
index 721f2d476..e643144ff 100644
--- a/Eigen/src/Core/MapBase.h
+++ b/Eigen/src/Core/MapBase.h
@@ -173,8 +173,14 @@ template<typename Derived> class MapBase
using Base::operator=;
using Base::operator*=;
- using Base::operator+=;
- using Base::operator-=;
+
+ template<typename Lhs,typename Rhs>
+ Derived& operator+=(const Flagged<Product<Lhs,Rhs,CacheFriendlyProduct>, 0, EvalBeforeNestingBit | EvalBeforeAssigningBit>& other)
+ { return Base::operator+=(other); }
+
+ template<typename Lhs,typename Rhs>
+ Derived& operator-=(const Flagged<Product<Lhs,Rhs,CacheFriendlyProduct>, 0, EvalBeforeNestingBit | EvalBeforeAssigningBit>& other)
+ { return Base::operator-=(other); }
template<typename OtherDerived>
Derived& operator+=(const MatrixBase<OtherDerived>& other)
diff --git a/Eigen/src/Core/Matrix.h b/Eigen/src/Core/Matrix.h
index d5c508128..f21364c70 100644
--- a/Eigen/src/Core/Matrix.h
+++ b/Eigen/src/Core/Matrix.h
@@ -571,10 +571,16 @@ class Matrix
template<typename OtherDerived>
EIGEN_STRONG_INLINE Matrix& _set(const MatrixBase<OtherDerived>& other)
{
- _resize_to_match(other);
- return Base::operator=(other);
+ _set_selector(other.derived(), typename ei_meta_if<(int(OtherDerived::Flags) & EvalBeforeAssigningBit), ei_meta_true, ei_meta_false>::ret());
+ return *this;
}
+ template<typename OtherDerived>
+ EIGEN_STRONG_INLINE void _set_selector(const OtherDerived& other, const ei_meta_true&) { _set_noalias(other.eval()); }
+
+ template<typename OtherDerived>
+ EIGEN_STRONG_INLINE void _set_selector(const OtherDerived& other, const ei_meta_false&) { _set_noalias(other); }
+
/** \internal Like _set() but additionally makes the assumption that no aliasing effect can happen (which
* is the case when creating a new matrix) so one can enforce lazy evaluation.
*
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index 3df518d1a..ac2b2532f 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -437,6 +437,7 @@ template<typename Derived> class MatrixBase
RealScalar norm() const;
RealScalar stableNorm() const;
RealScalar blueNorm() const;
+ RealScalar hypotNorm() const;
const PlainMatrixType normalized() const;
void normalize();
diff --git a/Eigen/src/Core/NumTraits.h b/Eigen/src/Core/NumTraits.h
index dec07a036..24afe54c5 100644
--- a/Eigen/src/Core/NumTraits.h
+++ b/Eigen/src/Core/NumTraits.h
@@ -70,9 +70,7 @@ template<> struct NumTraits<float>
HasFloatingPoint = 1,
ReadCost = 1,
AddCost = 1,
- MulCost = 1,
- Base = 2,
- Mantissa = 23
+ MulCost = 1
};
};
@@ -85,9 +83,7 @@ template<> struct NumTraits<double>
HasFloatingPoint = 1,
ReadCost = 1,
AddCost = 1,
- MulCost = 1,
- Base = 2,
- Mantissa = 52
+ MulCost = 1
};
};
diff --git a/Eigen/src/Core/Product.h b/Eigen/src/Core/Product.h
index 402f597d2..b46440ec0 100644
--- a/Eigen/src/Core/Product.h
+++ b/Eigen/src/Core/Product.h
@@ -287,6 +287,18 @@ MatrixBase<Derived>::operator*(const MatrixBase<OtherDerived> &other) const
return typename ProductReturnType<Derived,OtherDerived>::Type(derived(), other.derived());
}
+/** replaces \c *this by \c *this * \a other.
+ *
+ * \returns a reference to \c *this
+ */
+template<typename Derived>
+template<typename OtherDerived>
+inline Derived &
+MatrixBase<Derived>::operator*=(const MatrixBase<OtherDerived> &other)
+{
+ return derived() = derived() * other.derived();
+}
+
/***************************************************************************
* Normal product .coeff() implementation (with meta-unrolling)
***************************************************************************/
diff --git a/Eigen/src/Core/StableNorm.h b/Eigen/src/Core/StableNorm.h
new file mode 100644
index 000000000..22809633d
--- /dev/null
+++ b/Eigen/src/Core/StableNorm.h
@@ -0,0 +1,194 @@
+// This file is part of Eigen, a lightweight C++ template library
+// for linear algebra.
+//
+// Copyright (C) 2009 Gael Guennebaud <g.gael@free.fr>
+//
+// Eigen is free software; you can redistribute it and/or
+// modify it under the terms of the GNU Lesser General Public
+// License as published by the Free Software Foundation; either
+// version 3 of the License, or (at your option) any later version.
+//
+// Alternatively, you can redistribute it and/or
+// modify it under the terms of the GNU General Public License as
+// published by the Free Software Foundation; either version 2 of
+// the License, or (at your option) any later version.
+//
+// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
+// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
+// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
+// GNU General Public License for more details.
+//
+// You should have received a copy of the GNU Lesser General Public
+// License and a copy of the GNU General Public License along with
+// Eigen. If not, see <http://www.gnu.org/licenses/>.
+
+#ifndef EIGEN_STABLENORM_H
+#define EIGEN_STABLENORM_H
+
+template<typename ExpressionType, typename Scalar>
+inline void ei_stable_norm_kernel(const ExpressionType& bl, Scalar& ssq, Scalar& scale, Scalar& invScale)
+{
+ Scalar max = bl.cwise().abs().maxCoeff();
+ if (max>scale)
+ {
+ ssq = ssq * ei_abs2(scale/max);
+ scale = max;
+ invScale = Scalar(1)/scale;
+ }
+ // TODO if the max is much much smaller than the current scale,
+ // then we can neglect this sub vector
+ ssq += (bl*invScale).squaredNorm();
+}
+
+/** \returns the \em l2 norm of \c *this avoiding underflow and overflow.
+ * This version use a blockwise two passes algorithm:
+ * 1 - find the absolute largest coefficient \c s
+ * 2 - compute \f$ s \Vert \frac{*this}{s} \Vert \f$ in a standard way
+ *
+ * For architecture/scalar types supporting vectorization, this version
+ * is faster than blueNorm(). Otherwise the blueNorm() is much faster.
+ *
+ * \sa norm(), blueNorm(), hypotNorm()
+ */
+template<typename Derived>
+inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real
+MatrixBase<Derived>::stableNorm() const
+{
+ const int blockSize = 4096;
+ RealScalar scale = 0;
+ RealScalar invScale;
+ RealScalar ssq = 0; // sum of square
+ enum {
+ Alignment = (int(Flags)&DirectAccessBit) || (int(Flags)&AlignedBit) ? ForceAligned : AsRequested
+ };
+ int n = size();
+ int bi=0;
+ if ((int(Flags)&DirectAccessBit) && !(int(Flags)&AlignedBit))
+ {
+ bi = ei_alignmentOffset(&const_cast_derived().coeffRef(0), n);
+ if (bi>0)
+ ei_stable_norm_kernel(start(bi), ssq, scale, invScale);
+ }
+ for (; bi<n; bi+=blockSize)
+ ei_stable_norm_kernel(VectorBlock<Derived,Dynamic,Alignment>(derived(),bi,std::min(blockSize, n - bi)), ssq, scale, invScale);
+ return scale * ei_sqrt(ssq);
+}
+
+/** \returns the \em l2 norm of \c *this using the Blue's algorithm.
+ * A Portable Fortran Program to Find the Euclidean Norm of a Vector,
+ * ACM TOMS, Vol 4, Issue 1, 1978.
+ *
+ * For architecture/scalar types without vectorization, this version
+ * is much faster than stableNorm(). Otherwise the stableNorm() is faster.
+ *
+ * \sa norm(), stableNorm(), hypotNorm()
+ */
+template<typename Derived>
+inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real
+MatrixBase<Derived>::blueNorm() const
+{
+ static int nmax = -1;
+ static RealScalar b1, b2, s1m, s2m, overfl, rbig, relerr;
+ if(nmax <= 0)
+ {
+ int nbig, ibeta, it, iemin, iemax, iexp;
+ RealScalar abig, eps;
+ // This program calculates the machine-dependent constants
+ // bl, b2, slm, s2m, relerr overfl, nmax
+ // from the "basic" machine-dependent numbers
+ // nbig, ibeta, it, iemin, iemax, rbig.
+ // The following define the basic machine-dependent constants.
+ // For portability, the PORT subprograms "ilmaeh" and "rlmach"
+ // are used. For any specific computer, each of the assignment
+ // statements can be replaced
+ nbig = std::numeric_limits<int>::max(); // largest integer
+ ibeta = std::numeric_limits<RealScalar>::radix; // base for floating-point numbers
+ it = std::numeric_limits<RealScalar>::digits; // number of base-beta digits in mantissa
+ iemin = std::numeric_limits<RealScalar>::min_exponent; // minimum exponent
+ iemax = std::numeric_limits<RealScalar>::max_exponent; // maximum exponent
+ rbig = std::numeric_limits<RealScalar>::max(); // largest floating-point number
+
+ // Check the basic machine-dependent constants.
+ if(iemin > 1 - 2*it || 1+it>iemax || (it==2 && ibeta<5)
+ || (it<=4 && ibeta <= 3 ) || it<2)
+ {
+ ei_assert(false && "the algorithm cannot be guaranteed on this computer");
+ }
+ iexp = -((1-iemin)/2);
+ b1 = std::pow(double(ibeta),iexp); // lower boundary of midrange
+ iexp = (iemax + 1 - it)/2;
+ b2 = std::pow(double(ibeta),iexp); // upper boundary of midrange
+
+ iexp = (2-iemin)/2;
+ s1m = std::pow(double(ibeta),iexp); // scaling factor for lower range
+ iexp = - ((iemax+it)/2);
+ s2m = std::pow(double(ibeta),iexp); // scaling factor for upper range
+
+ overfl = rbig*s2m; // overfow boundary for abig
+ eps = std::pow(double(ibeta), 1-it);
+ relerr = ei_sqrt(eps); // tolerance for neglecting asml
+ abig = 1.0/eps - 1.0;
+ if (RealScalar(nbig)>abig) nmax = abig; // largest safe n
+ else nmax = nbig;
+ }
+ int n = size();
+ RealScalar ab2 = b2 / RealScalar(n);
+ RealScalar asml = RealScalar(0);
+ RealScalar amed = RealScalar(0);
+ RealScalar abig = RealScalar(0);
+ for(int j=0; j<n; ++j)
+ {
+ RealScalar ax = ei_abs(coeff(j));
+ if(ax > ab2) abig += ei_abs2(ax*s2m);
+ else if(ax < b1) asml += ei_abs2(ax*s1m);
+ else amed += ei_abs2(ax);
+ }
+ if(abig > RealScalar(0))
+ {
+ abig = ei_sqrt(abig);
+ if(abig > overfl)
+ {
+ ei_assert(false && "overflow");
+ return rbig;
+ }
+ if(amed > RealScalar(0))
+ {
+ abig = abig/s2m;
+ amed = ei_sqrt(amed);
+ }
+ else
+ return abig/s2m;
+ }
+ else if(asml > RealScalar(0))
+ {
+ if (amed > RealScalar(0))
+ {
+ abig = ei_sqrt(amed);
+ amed = ei_sqrt(asml) / s1m;
+ }
+ else
+ return ei_sqrt(asml)/s1m;
+ }
+ else
+ return ei_sqrt(amed);
+ asml = std::min(abig, amed);
+ abig = std::max(abig, amed);
+ if(asml <= abig*relerr)
+ return abig;
+ else
+ return abig * ei_sqrt(RealScalar(1) + ei_abs2(asml/abig));
+}
+
+/** \returns the \em l2 norm of \c *this avoiding undeflow and overflow.
+ * This version use a concatenation of hypot() calls, and it is very slow.
+ *
+ * \sa norm(), stableNorm()
+ */
+template<typename Derived>
+inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real
+MatrixBase<Derived>::hypotNorm() const
+{
+ return this->cwise().abs().redux(ei_scalar_hypot_op<RealScalar>());
+}
+
+#endif // EIGEN_STABLENORM_H