aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
-rw-r--r--Eigen/src/Core/Block.h22
-rw-r--r--Eigen/src/Core/Coeffs.h32
-rw-r--r--Eigen/src/Core/CommaInitializer.h14
-rw-r--r--Eigen/src/Core/CwiseBinaryOp.h2
-rw-r--r--Eigen/src/Core/DiagonalMatrix.h2
-rw-r--r--Eigen/src/Core/Dot.h2
-rw-r--r--Eigen/src/Core/Fuzzy.h4
-rw-r--r--Eigen/src/Core/Identity.h2
-rw-r--r--Eigen/src/Core/Map.h6
-rw-r--r--Eigen/src/Core/MathFunctions.h10
-rw-r--r--Eigen/src/Core/Matrix.h22
-rw-r--r--Eigen/src/Core/MatrixBase.h13
-rw-r--r--Eigen/src/Core/Minor.h2
-rw-r--r--Eigen/src/Core/Ones.h4
-rw-r--r--Eigen/src/Core/OperatorEquals.h6
-rw-r--r--Eigen/src/Core/Product.h6
-rw-r--r--Eigen/src/Core/Random.h4
-rw-r--r--Eigen/src/Core/Swap.h11
-rw-r--r--Eigen/src/Core/Util.h10
-rw-r--r--Eigen/src/Core/Visitor.h5
-rw-r--r--Eigen/src/Core/Zero.h4
-rw-r--r--doc/echelon.cpp20
22 files changed, 104 insertions, 99 deletions
diff --git a/Eigen/src/Core/Block.h b/Eigen/src/Core/Block.h
index c83cdc0de..6d39efe7e 100644
--- a/Eigen/src/Core/Block.h
+++ b/Eigen/src/Core/Block.h
@@ -90,7 +90,7 @@ template<typename MatrixType, int BlockRows, int BlockCols> class Block
m_blockRows(matrix.rows()), // if it is a row, then m_blockRows has a fixed-size of 1, so no pb to try to overwrite it
m_blockCols(matrix.cols()) // same for m_blockCols
{
- assert( (i>=0) && (
+ ei_assert( (i>=0) && (
((BlockRows==1) && (BlockCols==MatrixType::ColsAtCompileTime) && i<matrix.rows())
||((BlockRows==MatrixType::RowsAtCompileTime) && (BlockCols==1) && i<matrix.cols())));
}
@@ -100,8 +100,8 @@ template<typename MatrixType, int BlockRows, int BlockCols> class Block
Block(const MatrixType& matrix, int startRow, int startCol)
: m_matrix(matrix), m_startRow(startRow), m_startCol(startCol)
{
- assert(RowsAtCompileTime!=Dynamic && RowsAtCompileTime!=Dynamic);
- assert(startRow >= 0 && BlockRows >= 1 && startRow + BlockRows <= matrix.rows()
+ ei_assert(RowsAtCompileTime!=Dynamic && RowsAtCompileTime!=Dynamic);
+ ei_assert(startRow >= 0 && BlockRows >= 1 && startRow + BlockRows <= matrix.rows()
&& startCol >= 0 && BlockCols >= 1 && startCol + BlockCols <= matrix.cols());
}
@@ -113,9 +113,9 @@ template<typename MatrixType, int BlockRows, int BlockCols> class Block
: m_matrix(matrix), m_startRow(startRow), m_startCol(startCol),
m_blockRows(blockRows), m_blockCols(blockCols)
{
- assert((RowsAtCompileTime==Dynamic || RowsAtCompileTime==1)
+ ei_assert((RowsAtCompileTime==Dynamic || RowsAtCompileTime==1)
&& (ColsAtCompileTime==Dynamic || ColsAtCompileTime==1));
- assert(startRow >= 0 && blockRows >= 1 && startRow + blockRows <= matrix.rows()
+ ei_assert(startRow >= 0 && blockRows >= 1 && startRow + blockRows <= matrix.rows()
&& startCol >= 0 && blockCols >= 1 && startCol + blockCols <= matrix.cols());
}
@@ -197,7 +197,7 @@ template<typename Derived>
Block<Derived> MatrixBase<Derived>
::block(int start, int size)
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
return Block<Derived>(derived(), RowsAtCompileTime == 1 ? 0 : start,
ColsAtCompileTime == 1 ? 0 : start,
RowsAtCompileTime == 1 ? 1 : size,
@@ -209,7 +209,7 @@ template<typename Derived>
const Block<Derived> MatrixBase<Derived>
::block(int start, int size) const
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
return Block<Derived>(derived(), RowsAtCompileTime == 1 ? 0 : start,
ColsAtCompileTime == 1 ? 0 : start,
RowsAtCompileTime == 1 ? 1 : size,
@@ -235,7 +235,7 @@ template<typename Derived>
Block<Derived> MatrixBase<Derived>
::start(int size)
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
return Block<Derived>(derived(), 0, 0,
RowsAtCompileTime == 1 ? 1 : size,
ColsAtCompileTime == 1 ? 1 : size);
@@ -246,7 +246,7 @@ template<typename Derived>
const Block<Derived> MatrixBase<Derived>
::start(int size) const
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
return Block<Derived>(derived(), 0, 0,
RowsAtCompileTime == 1 ? 1 : size,
ColsAtCompileTime == 1 ? 1 : size);
@@ -271,7 +271,7 @@ template<typename Derived>
Block<Derived> MatrixBase<Derived>
::end(int size)
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
return Block<Derived>(derived(),
RowsAtCompileTime == 1 ? 0 : rows() - size,
ColsAtCompileTime == 1 ? 0 : cols() - size,
@@ -284,7 +284,7 @@ template<typename Derived>
const Block<Derived> MatrixBase<Derived>
::end(int size) const
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
return Block<Derived>(derived(),
RowsAtCompileTime == 1 ? 0 : rows() - size,
ColsAtCompileTime == 1 ? 0 : cols() - size,
diff --git a/Eigen/src/Core/Coeffs.h b/Eigen/src/Core/Coeffs.h
index 905093f7b..d0b211314 100644
--- a/Eigen/src/Core/Coeffs.h
+++ b/Eigen/src/Core/Coeffs.h
@@ -43,7 +43,7 @@ template<typename Derived>
typename ei_traits<Derived>::Scalar MatrixBase<Derived>
::coeff(int row, int col) const
{
- eigen_internal_assert(row >= 0 && row < rows()
+ ei_internal_assert(row >= 0 && row < rows()
&& col >= 0 && col < cols());
return derived()._coeff(row, col);
}
@@ -56,7 +56,7 @@ template<typename Derived>
typename ei_traits<Derived>::Scalar MatrixBase<Derived>
::operator()(int row, int col) const
{
- assert(row >= 0 && row < rows()
+ ei_assert(row >= 0 && row < rows()
&& col >= 0 && col < cols());
return derived()._coeff(row, col);
}
@@ -79,7 +79,7 @@ template<typename Derived>
typename ei_traits<Derived>::Scalar& MatrixBase<Derived>
::coeffRef(int row, int col)
{
- eigen_internal_assert(row >= 0 && row < rows()
+ ei_internal_assert(row >= 0 && row < rows()
&& col >= 0 && col < cols());
return derived()._coeffRef(row, col);
}
@@ -92,7 +92,7 @@ template<typename Derived>
typename ei_traits<Derived>::Scalar& MatrixBase<Derived>
::operator()(int row, int col)
{
- assert(row >= 0 && row < rows()
+ ei_assert(row >= 0 && row < rows()
&& col >= 0 && col < cols());
return derived()._coeffRef(row, col);
}
@@ -115,15 +115,15 @@ template<typename Derived>
typename ei_traits<Derived>::Scalar MatrixBase<Derived>
::coeff(int index) const
{
- eigen_internal_assert(IsVectorAtCompileTime);
+ ei_internal_assert(IsVectorAtCompileTime);
if(RowsAtCompileTime == 1)
{
- eigen_internal_assert(index >= 0 && index < cols());
+ ei_internal_assert(index >= 0 && index < cols());
return coeff(0, index);
}
else
{
- eigen_internal_assert(index >= 0 && index < rows());
+ ei_internal_assert(index >= 0 && index < rows());
return coeff(index, 0);
}
}
@@ -139,15 +139,15 @@ template<typename Derived>
typename ei_traits<Derived>::Scalar MatrixBase<Derived>
::operator[](int index) const
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
if(RowsAtCompileTime == 1)
{
- assert(index >= 0 && index < cols());
+ ei_assert(index >= 0 && index < cols());
return coeff(0, index);
}
else
{
- assert(index >= 0 && index < rows());
+ ei_assert(index >= 0 && index < rows());
return coeff(index, 0);
}
}
@@ -170,15 +170,15 @@ template<typename Derived>
typename ei_traits<Derived>::Scalar& MatrixBase<Derived>
::coeffRef(int index)
{
- eigen_internal_assert(IsVectorAtCompileTime);
+ ei_internal_assert(IsVectorAtCompileTime);
if(RowsAtCompileTime == 1)
{
- eigen_internal_assert(index >= 0 && index < cols());
+ ei_internal_assert(index >= 0 && index < cols());
return coeffRef(0, index);
}
else
{
- eigen_internal_assert(index >= 0 && index < rows());
+ ei_internal_assert(index >= 0 && index < rows());
return coeffRef(index, 0);
}
}
@@ -193,15 +193,15 @@ template<typename Derived>
typename ei_traits<Derived>::Scalar& MatrixBase<Derived>
::operator[](int index)
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
if(RowsAtCompileTime == 1)
{
- assert(index >= 0 && index < cols());
+ ei_assert(index >= 0 && index < cols());
return coeffRef(0, index);
}
else
{
- assert(index >= 0 && index < rows());
+ ei_assert(index >= 0 && index < rows());
return coeffRef(index, 0);
}
}
diff --git a/Eigen/src/Core/CommaInitializer.h b/Eigen/src/Core/CommaInitializer.h
index c7a292437..656ee3002 100644
--- a/Eigen/src/Core/CommaInitializer.h
+++ b/Eigen/src/Core/CommaInitializer.h
@@ -52,12 +52,12 @@ struct MatrixBase<Derived>::CommaInitializer
m_row+=m_currentBlockRows;
m_col = 0;
m_currentBlockRows = 1;
- assert(m_row<m_matrix.rows()
+ ei_assert(m_row<m_matrix.rows()
&& "Too many rows passed to MatrixBase::operator<<");
}
- assert(m_col<m_matrix.cols()
+ ei_assert(m_col<m_matrix.cols()
&& "Too many coefficients passed to MatrixBase::operator<<");
- assert(m_currentBlockRows==1);
+ ei_assert(m_currentBlockRows==1);
m_matrix.coeffRef(m_row, m_col++) = s;
return *this;
}
@@ -70,12 +70,12 @@ struct MatrixBase<Derived>::CommaInitializer
m_row+=m_currentBlockRows;
m_col = 0;
m_currentBlockRows = other.rows();
- assert(m_row+m_currentBlockRows<=m_matrix.rows()
+ ei_assert(m_row+m_currentBlockRows<=m_matrix.rows()
&& "Too many rows passed to MatrixBase::operator<<");
}
- assert(m_col<m_matrix.cols()
+ ei_assert(m_col<m_matrix.cols()
&& "Too many coefficients passed to MatrixBase::operator<<");
- assert(m_currentBlockRows==other.rows());
+ ei_assert(m_currentBlockRows==other.rows());
if (OtherDerived::RowsAtCompileTime>0 && OtherDerived::ColsAtCompileTime>0)
m_matrix.block< (OtherDerived::RowsAtCompileTime>0?OtherDerived::RowsAtCompileTime:1) ,
(OtherDerived::ColsAtCompileTime>0?OtherDerived::ColsAtCompileTime:1) >(m_row, m_col) = other;
@@ -87,7 +87,7 @@ struct MatrixBase<Derived>::CommaInitializer
~CommaInitializer(void)
{
- assert((m_row+m_currentBlockRows) == m_matrix.rows()
+ ei_assert((m_row+m_currentBlockRows) == m_matrix.rows()
&& m_col == m_matrix.cols()
&& "Too few coefficients passed to Matrix::operator<<");
}
diff --git a/Eigen/src/Core/CwiseBinaryOp.h b/Eigen/src/Core/CwiseBinaryOp.h
index 1ab759241..ff7e8661c 100644
--- a/Eigen/src/Core/CwiseBinaryOp.h
+++ b/Eigen/src/Core/CwiseBinaryOp.h
@@ -74,7 +74,7 @@ class CwiseBinaryOp : ei_no_assignment_operator,
CwiseBinaryOp(const Lhs& lhs, const Rhs& rhs, const BinaryOp& func = BinaryOp())
: m_lhs(lhs), m_rhs(rhs), m_functor(func)
{
- assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
+ ei_assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
}
private:
diff --git a/Eigen/src/Core/DiagonalMatrix.h b/Eigen/src/Core/DiagonalMatrix.h
index 79fa562d9..36d2c656b 100644
--- a/Eigen/src/Core/DiagonalMatrix.h
+++ b/Eigen/src/Core/DiagonalMatrix.h
@@ -60,7 +60,7 @@ class DiagonalMatrix : ei_no_assignment_operator,
DiagonalMatrix(const CoeffsVectorType& coeffs) : m_coeffs(coeffs)
{
- assert(CoeffsVectorType::IsVectorAtCompileTime
+ ei_assert(CoeffsVectorType::IsVectorAtCompileTime
&& coeffs.size() > 0);
}
diff --git a/Eigen/src/Core/Dot.h b/Eigen/src/Core/Dot.h
index 33f69f987..248fbca79 100644
--- a/Eigen/src/Core/Dot.h
+++ b/Eigen/src/Core/Dot.h
@@ -72,7 +72,7 @@ template<typename OtherDerived>
typename ei_traits<Derived>::Scalar
MatrixBase<Derived>::dot(const MatrixBase<OtherDerived>& other) const
{
- assert(IsVectorAtCompileTime
+ ei_assert(IsVectorAtCompileTime
&& OtherDerived::IsVectorAtCompileTime
&& size() == other.size());
Scalar res;
diff --git a/Eigen/src/Core/Fuzzy.h b/Eigen/src/Core/Fuzzy.h
index 0f601a8f5..4ebe9d2a7 100644
--- a/Eigen/src/Core/Fuzzy.h
+++ b/Eigen/src/Core/Fuzzy.h
@@ -48,7 +48,7 @@ bool MatrixBase<Derived>::isApprox(
typename NumTraits<Scalar>::Real prec
) const
{
- assert(rows() == other.rows() && cols() == other.cols());
+ ei_assert(rows() == other.rows() && cols() == other.cols());
if(IsVectorAtCompileTime)
{
return((*this - other).norm2() <= std::min(norm2(), other.norm2()) * prec * prec);
@@ -109,7 +109,7 @@ bool MatrixBase<Derived>::isMuchSmallerThan(
typename NumTraits<Scalar>::Real prec
) const
{
- assert(rows() == other.rows() && cols() == other.cols());
+ ei_assert(rows() == other.rows() && cols() == other.cols());
if(IsVectorAtCompileTime)
{
return(norm2() <= other.norm2() * prec * prec);
diff --git a/Eigen/src/Core/Identity.h b/Eigen/src/Core/Identity.h
index 129786e4b..0b39a4274 100644
--- a/Eigen/src/Core/Identity.h
+++ b/Eigen/src/Core/Identity.h
@@ -52,7 +52,7 @@ template<typename MatrixType> class Identity : ei_no_assignment_operator,
Identity(int rows, int cols) : m_rows(rows), m_cols(cols)
{
- assert(rows > 0
+ ei_assert(rows > 0
&& (RowsAtCompileTime == Dynamic || RowsAtCompileTime == rows)
&& cols > 0
&& (ColsAtCompileTime == Dynamic || ColsAtCompileTime == cols));
diff --git a/Eigen/src/Core/Map.h b/Eigen/src/Core/Map.h
index d0f9cd940..5d08ef542 100644
--- a/Eigen/src/Core/Map.h
+++ b/Eigen/src/Core/Map.h
@@ -81,7 +81,7 @@ template<typename MatrixType> class Map
public:
Map(const Scalar* data, int rows, int cols) : m_data(data), m_rows(rows), m_cols(cols)
{
- assert(rows > 0
+ ei_assert(rows > 0
&& (RowsAtCompileTime == Dynamic || RowsAtCompileTime == rows)
&& cols > 0
&& (ColsAtCompileTime == Dynamic || ColsAtCompileTime == cols));
@@ -107,7 +107,7 @@ template<typename _Scalar, int _Rows, int _Cols, int _StorageOrder, int _MaxRows
const Map<Matrix<_Scalar, _Rows, _Cols, _StorageOrder, _MaxRows, _MaxCols> >
Matrix<_Scalar, _Rows, _Cols, _StorageOrder, _MaxRows, _MaxCols>::map(const Scalar* data, int size)
{
- assert(_Cols == 1 || _Rows ==1);
+ ei_assert(_Cols == 1 || _Rows ==1);
if(_Cols == 1)
return Map<Matrix>(data, size, 1);
else
@@ -156,7 +156,7 @@ template<typename _Scalar, int _Rows, int _Cols, int _StorageOrder, int _MaxRows
Map<Matrix<_Scalar, _Rows, _Cols, _StorageOrder, _MaxRows, _MaxCols> >
Matrix<_Scalar, _Rows, _Cols, _StorageOrder, _MaxRows, _MaxCols>::map(Scalar* data, int size)
{
- assert(_Cols == 1 || _Rows ==1);
+ ei_assert(_Cols == 1 || _Rows ==1);
if(_Cols == 1)
return Map<Matrix>(data, size, 1);
else
diff --git a/Eigen/src/Core/MathFunctions.h b/Eigen/src/Core/MathFunctions.h
index a1a54cc38..df30c5ad4 100644
--- a/Eigen/src/Core/MathFunctions.h
+++ b/Eigen/src/Core/MathFunctions.h
@@ -35,11 +35,11 @@ inline int ei_imag(int) { return 0; }
inline int ei_conj(int x) { return x; }
inline int ei_abs(int x) { return abs(x); }
inline int ei_abs2(int x) { return x*x; }
-inline int ei_sqrt(int) { assert(false); return 0; }
-inline int ei_exp(int) { assert(false); return 0; }
-inline int ei_log(int) { assert(false); return 0; }
-inline int ei_sin(int) { assert(false); return 0; }
-inline int ei_cos(int) { assert(false); return 0; }
+inline int ei_sqrt(int) { ei_assert(false); return 0; }
+inline int ei_exp(int) { ei_assert(false); return 0; }
+inline int ei_log(int) { ei_assert(false); return 0; }
+inline int ei_sin(int) { ei_assert(false); return 0; }
+inline int ei_cos(int) { ei_assert(false); return 0; }
#if (defined __ICC) || (defined __GNUC__ && (__GNUC__<4 || __GNUC_MINOR__<3))
inline int ei_pow(int x, int y) { return int(std::pow(double(x), y)); }
#else
diff --git a/Eigen/src/Core/Matrix.h b/Eigen/src/Core/Matrix.h
index 9452906d0..60365f625 100644
--- a/Eigen/src/Core/Matrix.h
+++ b/Eigen/src/Core/Matrix.h
@@ -131,7 +131,7 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols,
void resize(int rows, int cols)
{
- assert(rows > 0
+ ei_assert(rows > 0
&& (MaxRowsAtCompileTime == Dynamic || MaxRowsAtCompileTime >= rows)
&& (RowsAtCompileTime == Dynamic || RowsAtCompileTime == rows)
&& cols > 0
@@ -153,12 +153,12 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols,
{
if(RowsAtCompileTime == 1)
{
- assert(other.isVector());
+ ei_assert(other.isVector());
resize(1, other.size());
}
else if(ColsAtCompileTime == 1)
{
- assert(other.isVector());
+ ei_assert(other.isVector());
resize(other.size(), 1);
}
else resize(other.rows(), other.cols());
@@ -191,7 +191,7 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols,
*/
explicit Matrix()
{
- assert(RowsAtCompileTime > 0 && ColsAtCompileTime > 0);
+ ei_assert(RowsAtCompileTime > 0 && ColsAtCompileTime > 0);
}
/** Constructs a vector or row-vector with given dimension. \only_for_vectors
@@ -202,8 +202,8 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols,
*/
explicit Matrix(int dim) : m_storage(dim, RowsAtCompileTime == 1 ? 1 : dim, ColsAtCompileTime == 1 ? 1 : dim)
{
- assert(dim > 0);
- assert((RowsAtCompileTime == 1
+ ei_assert(dim > 0);
+ ei_assert((RowsAtCompileTime == 1
&& (ColsAtCompileTime == Dynamic || ColsAtCompileTime == dim))
|| (ColsAtCompileTime == 1
&& (RowsAtCompileTime == Dynamic || RowsAtCompileTime == dim)));
@@ -229,14 +229,14 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols,
}
else
{
- assert(x > 0 && (RowsAtCompileTime == Dynamic || RowsAtCompileTime == x)
+ ei_assert(x > 0 && (RowsAtCompileTime == Dynamic || RowsAtCompileTime == x)
&& y > 0 && (ColsAtCompileTime == Dynamic || ColsAtCompileTime == y));
}
}
/** constructs an initialized 2D vector with given coefficients */
Matrix(const float& x, const float& y)
{
- assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 2)
+ ei_assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 2)
|| (RowsAtCompileTime == 2 && ColsAtCompileTime == 1));
m_storage.data()[0] = x;
m_storage.data()[1] = y;
@@ -244,7 +244,7 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols,
/** constructs an initialized 2D vector with given coefficients */
Matrix(const double& x, const double& y)
{
- assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 2)
+ ei_assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 2)
|| (RowsAtCompileTime == 2 && ColsAtCompileTime == 1));
m_storage.data()[0] = x;
m_storage.data()[1] = y;
@@ -252,7 +252,7 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols,
/** constructs an initialized 3D vector with given coefficients */
Matrix(const Scalar& x, const Scalar& y, const Scalar& z)
{
- assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 3)
+ ei_assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 3)
|| (RowsAtCompileTime == 3 && ColsAtCompileTime == 1));
m_storage.data()[0] = x;
m_storage.data()[1] = y;
@@ -261,7 +261,7 @@ class Matrix : public MatrixBase<Matrix<_Scalar, _Rows, _Cols,
/** constructs an initialized 4D vector with given coefficients */
Matrix(const Scalar& x, const Scalar& y, const Scalar& z, const Scalar& w)
{
- assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 4)
+ ei_assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 4)
|| (RowsAtCompileTime == 4 && ColsAtCompileTime == 1));
m_storage.data()[0] = x;
m_storage.data()[1] = y;
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index d2f415741..fa3a09dd6 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -343,13 +343,6 @@ template<typename Derived> class MatrixBase
RealScalar prec = precision<Scalar>()) const;
bool isOrtho(RealScalar prec = precision<Scalar>()) const;
- /** puts in *row and *col the location of the coefficient of *this
- * which has the biggest absolute value.
- */
- void findBiggestCoeff(int *row, int *col) const
- { (*this).cwiseAbs().maxCoeff(row, col); }
- //@}
-
/// \name Special functions
//@{
template<typename NewType>
@@ -358,12 +351,6 @@ template<typename Derived> class MatrixBase
const Eval<Derived> eval() const EIGEN_ALWAYS_INLINE;
const EvalOMP<Derived> evalOMP() const EIGEN_ALWAYS_INLINE;
- /** swaps *this with the expression \a other.
- *
- * \note \a other is only marked const because I couln't find another way
- * to get g++ 4.2 to accept that template parameter resolution. It gets const_cast'd
- * of course. TODO: get rid of const here.
- */
template<typename OtherDerived>
void swap(const MatrixBase<OtherDerived>& other);
//@}
diff --git a/Eigen/src/Core/Minor.h b/Eigen/src/Core/Minor.h
index 6c755bb74..3897c8bac 100644
--- a/Eigen/src/Core/Minor.h
+++ b/Eigen/src/Core/Minor.h
@@ -64,7 +64,7 @@ template<typename MatrixType> class Minor
int row, int col)
: m_matrix(matrix), m_row(row), m_col(col)
{
- assert(row >= 0 && row < matrix.rows()
+ ei_assert(row >= 0 && row < matrix.rows()
&& col >= 0 && col < matrix.cols());
}
diff --git a/Eigen/src/Core/Ones.h b/Eigen/src/Core/Ones.h
index 15ecc04e0..0d0f002bd 100644
--- a/Eigen/src/Core/Ones.h
+++ b/Eigen/src/Core/Ones.h
@@ -64,7 +64,7 @@ template<typename MatrixType> class Ones : ei_no_assignment_operator,
public:
Ones(int rows, int cols) : m_rows(rows), m_cols(cols)
{
- assert(rows > 0
+ ei_assert(rows > 0
&& (RowsAtCompileTime == Dynamic || RowsAtCompileTime == rows)
&& cols > 0
&& (ColsAtCompileTime == Dynamic || ColsAtCompileTime == cols));
@@ -114,7 +114,7 @@ const Ones<Derived> MatrixBase<Derived>::ones(int rows, int cols)
template<typename Derived>
const Ones<Derived> MatrixBase<Derived>::ones(int size)
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
if(RowsAtCompileTime == 1) return Ones<Derived>(1, size);
else return Ones<Derived>(size, 1);
}
diff --git a/Eigen/src/Core/OperatorEquals.h b/Eigen/src/Core/OperatorEquals.h
index 82ef13ead..babd51776 100644
--- a/Eigen/src/Core/OperatorEquals.h
+++ b/Eigen/src/Core/OperatorEquals.h
@@ -105,7 +105,7 @@ Derived& MatrixBase<Derived>
if(IsVectorAtCompileTime && OtherDerived::IsVectorAtCompileTime)
// copying a vector expression into a vector
{
- assert(size() == other.size());
+ ei_assert(size() == other.size());
if(EIGEN_UNROLLED_LOOPS
&& SizeAtCompileTime != Dynamic
&& SizeAtCompileTime <= EIGEN_UNROLLING_LIMIT)
@@ -120,7 +120,7 @@ Derived& MatrixBase<Derived>
}
else // copying a matrix expression into a matrix
{
- assert(rows() == other.rows() && cols() == other.cols());
+ ei_assert(rows() == other.rows() && cols() == other.cols());
if(EIGEN_UNROLLED_LOOPS
&& SizeAtCompileTime != Dynamic
&& SizeAtCompileTime <= EIGEN_UNROLLING_LIMIT)
@@ -148,7 +148,7 @@ Derived& MatrixBase<Derived>
coeffRef(i, j) = other.coeff(i, j);
}
}
- return (*this).derived();
+ return derived();
}
}
diff --git a/Eigen/src/Core/Product.h b/Eigen/src/Core/Product.h
index 28602d8a0..13a995ad0 100644
--- a/Eigen/src/Core/Product.h
+++ b/Eigen/src/Core/Product.h
@@ -106,7 +106,7 @@ template<typename Lhs, typename Rhs, int EvalMode> class Product : ei_no_assignm
Product(const Lhs& lhs, const Rhs& rhs)
: m_lhs(lhs), m_rhs(rhs)
{
- assert(lhs.cols() == rhs.rows());
+ ei_assert(lhs.cols() == rhs.rows());
}
/** \internal */
@@ -172,7 +172,7 @@ template<typename OtherDerived>
const Eval<Product<Derived, OtherDerived> >
MatrixBase<Derived>::operator*(const MatrixBase<OtherDerived> &other) const
{
- return (*this).lazyProduct(other).eval();
+ return lazyProduct(other).eval();
}
/** replaces \c *this by \c *this * \a other.
@@ -192,7 +192,7 @@ template<typename Derived1, typename Derived2>
Derived& MatrixBase<Derived>::operator=(const Product<Derived1,Derived2,CacheOptimal>& product)
{
product._cacheOptimalEval(*this);
- return (*this).derived();
+ return derived();
}
template<typename Lhs, typename Rhs, int EvalMode>
diff --git a/Eigen/src/Core/Random.h b/Eigen/src/Core/Random.h
index 188f7513a..f756c9649 100644
--- a/Eigen/src/Core/Random.h
+++ b/Eigen/src/Core/Random.h
@@ -65,7 +65,7 @@ template<typename MatrixType> class Random : ei_no_assignment_operator,
Random(int rows, int cols) : m_rows(rows), m_cols(cols)
{
- assert(rows > 0
+ ei_assert(rows > 0
&& (RowsAtCompileTime == Dynamic || RowsAtCompileTime == rows)
&& cols > 0
&& (ColsAtCompileTime == Dynamic || ColsAtCompileTime == cols));
@@ -117,7 +117,7 @@ template<typename Derived>
const Eval<Random<Derived> >
MatrixBase<Derived>::random(int size)
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
if(RowsAtCompileTime == 1) return Random<Derived>(1, size).eval();
else return Random<Derived>(size, 1).eval();
}
diff --git a/Eigen/src/Core/Swap.h b/Eigen/src/Core/Swap.h
index 5bd16530e..db096e6f2 100644
--- a/Eigen/src/Core/Swap.h
+++ b/Eigen/src/Core/Swap.h
@@ -25,6 +25,15 @@
#ifndef EIGEN_SWAP_H
#define EIGEN_SWAP_H
+/** swaps *this with the expression \a other.
+ *
+ * \note \a other is only marked const because I couln't find another way
+ * to get g++ (4.2 and 4.3) to accept that template parameter resolution.
+ * The problem seems to be that when swapping expressions as in
+ * m.row(i).swap(m.row(j)); the Row object returned by row(j) is a temporary
+ * and g++ doesn't dare to pass it by non-constant reference.
+ * It gets const_cast'd of course. TODO: get rid of const here.
+ */
template<typename Derived>
template<typename OtherDerived>
void MatrixBase<Derived>::swap(const MatrixBase<OtherDerived>& other)
@@ -35,7 +44,7 @@ void MatrixBase<Derived>::swap(const MatrixBase<OtherDerived>& other)
Scalar tmp;
if(IsVectorAtCompileTime)
{
- assert(OtherDerived::IsVectorAtCompileTime && size() == _other->size());
+ ei_assert(OtherDerived::IsVectorAtCompileTime && size() == _other->size());
for(int i = 0; i < size(); i++)
{
tmp = coeff(i);
diff --git a/Eigen/src/Core/Util.h b/Eigen/src/Core/Util.h
index 454edd4cb..b458fe4ae 100644
--- a/Eigen/src/Core/Util.h
+++ b/Eigen/src/Core/Util.h
@@ -49,10 +49,16 @@ EIGEN_USING_MATRIX_TYPEDEFS \
using Eigen::Matrix; \
using Eigen::MatrixBase;
+#ifdef EIGEN_NDEBUG
+#define ei_assert(x)
+#else
+#define ei_assert(x) assert(x)
+#endif
+
#ifdef EIGEN_INTERNAL_DEBUGGING
-#define eigen_internal_assert(x) assert(x);
+#define ei_internal_assert(x) ei_assert(x);
#else
-#define eigen_internal_assert(x)
+#define ei_internal_assert(x)
#endif
#ifdef NDEBUG
diff --git a/Eigen/src/Core/Visitor.h b/Eigen/src/Core/Visitor.h
index d5b6fb167..7e3a7b819 100644
--- a/Eigen/src/Core/Visitor.h
+++ b/Eigen/src/Core/Visitor.h
@@ -82,13 +82,12 @@ void MatrixBase<Derived>::visit(Visitor& visitor) const
SizeAtCompileTime : Dynamic>::run(derived(), visitor);
else
{
- Scalar res;
visitor.init(coeff(0,0), 0, 0);
for(int i = 1; i < rows(); i++)
- visitor(res, coeff(i, 0), i, 0);
+ visitor(coeff(i, 0), i, 0);
for(int j = 1; j < cols(); j++)
for(int i = 0; i < rows(); i++)
- visitor(res, coeff(i, j), i, j);
+ visitor(coeff(i, j), i, j);
}
}
diff --git a/Eigen/src/Core/Zero.h b/Eigen/src/Core/Zero.h
index 97bb8b9e9..33792f0b1 100644
--- a/Eigen/src/Core/Zero.h
+++ b/Eigen/src/Core/Zero.h
@@ -65,7 +65,7 @@ template<typename MatrixType> class Zero : ei_no_assignment_operator,
Zero(int rows, int cols) : m_rows(rows), m_cols(cols)
{
- assert(rows > 0
+ ei_assert(rows > 0
&& (RowsAtCompileTime == Dynamic || RowsAtCompileTime == rows)
&& cols > 0
&& (ColsAtCompileTime == Dynamic || ColsAtCompileTime == cols));
@@ -115,7 +115,7 @@ const Zero<Derived> MatrixBase<Derived>::zero(int rows, int cols)
template<typename Derived>
const Zero<Derived> MatrixBase<Derived>::zero(int size)
{
- assert(IsVectorAtCompileTime);
+ ei_assert(IsVectorAtCompileTime);
if(RowsAtCompileTime == 1) return Zero<Derived>(1, size);
else return Zero<Derived>(size, 1);
}
diff --git a/doc/echelon.cpp b/doc/echelon.cpp
index 5b07db421..cf15fda49 100644
--- a/doc/echelon.cpp
+++ b/doc/echelon.cpp
@@ -7,19 +7,22 @@ namespace Eigen {
template<typename Derived>
void echelon(MatrixBase<Derived>& m)
{
- const int N = std::min(m.rows(), m.cols());
-
- for(int k = 0; k < N; k++)
+ for(int k = 0; k < m.diagonal().size(); k++)
{
int rowOfBiggest, colOfBiggest;
- int cornerRows = m.rows()-k;
- int cornerCols = m.cols()-k;
+ int cornerRows = m.rows()-k, cornerCols = m.cols()-k;
m.corner(BottomRight, cornerRows, cornerCols)
- .findBiggestCoeff(&rowOfBiggest, &colOfBiggest);
+ .cwiseAbs()
+ .maxCoeff(&rowOfBiggest, &colOfBiggest);
m.row(k).swap(m.row(k+rowOfBiggest));
m.col(k).swap(m.col(k+colOfBiggest));
- for(int r = k+1; r < m.rows(); r++)
- m.row(r).end(cornerCols) -= m.row(k).end(cornerCols) * m(r,k) / m(k,k);
+ // important performance tip:
+ // in a complex expression such as below it can be very important to fine-tune
+ // exactly where evaluation occurs. The parentheses and .eval() below ensure
+ // that the quotient is computed only once, and that the evaluation caused
+ // by operator* occurs last.
+ m.corner(BottomRight, cornerRows-1, cornerCols)
+ -= m.col(k).end(cornerRows-1) * (m.row(k).end(cornerCols) / m(k,k)).eval();
}
}
@@ -65,6 +68,7 @@ int main(int, char **)
cout << "Here's the matrix m:" << endl << m << endl;
cout << "Now let's echelon m:" << endl;
+ for(int i = 0; i < 1000000; i++)
echelon(m);
cout << "Now m is:" << endl << m << endl;