aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2008-06-06 13:10:00 +0000
committerGravatar Benoit Jacob <jacob.benoit.1@gmail.com>2008-06-06 13:10:00 +0000
commit869394ee8b0faa11f03412a77035fbf6c5c66f95 (patch)
tree4525a93fdc237bc07306178d2751fb3b9da2c9b5
parent2126baf9dcdbfb27579522f42228316c4c309539 (diff)
fix some compile errors with gcc 4.3, some warnings, some documentation
-rw-r--r--Eigen/src/Core/MatrixBase.h6
-rw-r--r--Eigen/src/QR/SelfAdjointEigenSolver.h9
-rwxr-xr-xEigen/src/QR/Tridiagonalization.h4
-rw-r--r--test/sizeof.cpp2
4 files changed, 9 insertions, 12 deletions
diff --git a/Eigen/src/Core/MatrixBase.h b/Eigen/src/Core/MatrixBase.h
index ea693b7a9..7e22caa75 100644
--- a/Eigen/src/Core/MatrixBase.h
+++ b/Eigen/src/Core/MatrixBase.h
@@ -141,9 +141,6 @@ template<typename Derived> class MatrixBase : public ArrayBase<Derived>
* (e.g. int, float or double) then \a RealScalar is just the same as \a Scalar. If
* \a Scalar is \a std::complex<T> then RealScalar is \a T.
*
- * In fact, \a RealScalar is defined as follows:
- * \code typedef typename NumTraits<Scalar>::Real RealScalar; \endcode
- *
* \sa class NumTraits
*/
typedef typename NumTraits<Scalar>::Real RealScalar;
@@ -161,7 +158,6 @@ template<typename Derived> class MatrixBase : public ArrayBase<Derived>
* \sa rows(), cols(), IsVectorAtCompileTime. */
inline bool isVector() const { return rows()==1 || cols()==1; }
-
/** Represents a constant matrix */
typedef CwiseNullaryOp<ei_scalar_constant_op<Scalar>,Derived> ConstantReturnType;
/** Represents a vector block of a matrix */
@@ -433,7 +429,7 @@ template<typename Derived> class MatrixBase : public ArrayBase<Derived>
/** \returns number of elements to skip to pass from one row (resp. column) to another
* for a row-major (resp. column-major) matrix.
- * Combined with coeffRef() and the compile times flags, it allows a direct access to the data
+ * Combined with coeffRef() and the \ref flags flags, it allows a direct access to the data
* of the underlying matrix.
*/
inline int stride(void) const { return derived()._stride(); }
diff --git a/Eigen/src/QR/SelfAdjointEigenSolver.h b/Eigen/src/QR/SelfAdjointEigenSolver.h
index 227d71528..fd73714dc 100644
--- a/Eigen/src/QR/SelfAdjointEigenSolver.h
+++ b/Eigen/src/QR/SelfAdjointEigenSolver.h
@@ -47,7 +47,7 @@ template<typename _MatrixType> class SelfAdjointEigenSolver
typedef std::complex<RealScalar> Complex;
typedef Matrix<RealScalar, MatrixType::ColsAtCompileTime, 1> RealVectorType;
typedef Matrix<RealScalar, Dynamic, 1> RealVectorTypeX;
- typedef Tridiagonalization<MatrixType> Tridiagonalization;
+ typedef Tridiagonalization<MatrixType> TridiagonalizationType;
SelfAdjointEigenSolver(const MatrixType& matrix, bool computeEigenvectors = true)
: m_eivec(matrix.rows(), matrix.cols()),
@@ -124,8 +124,8 @@ void SelfAdjointEigenSolver<MatrixType>::compute(const MatrixType& matrix, bool
// the latter avoids multiple memory allocation when the same SelfAdjointEigenSolver is used multiple times...
// (same for diag and subdiag)
RealVectorType& diag = m_eivalues;
- typename Tridiagonalization::SubDiagonalType subdiag(n-1);
- Tridiagonalization::decomposeInPlace(m_eivec, diag, subdiag, computeEigenvectors);
+ typename TridiagonalizationType::SubDiagonalType subdiag(n-1);
+ TridiagonalizationType::decomposeInPlace(m_eivec, diag, subdiag, computeEigenvectors);
int end = n-1;
int start = 0;
@@ -191,10 +191,11 @@ template<typename Derived> struct ei_matrixNorm_selector<Derived, false>
static inline typename NumTraits<typename ei_traits<Derived>::Scalar>::Real
matrixNorm(const MatrixBase<Derived>& m)
{
+ typename Derived::Eval m_eval(m);
// FIXME if it is really guaranteed that the eigenvalues are already sorted,
// then we don't need to compute a maxCoeff() here, comparing the 1st and last ones is enough.
return ei_sqrt(
- (m*m.adjoint())
+ (m_eval*m_eval.adjoint())
.template marked<SelfAdjoint>()
.eigenvalues()
.maxCoeff()
diff --git a/Eigen/src/QR/Tridiagonalization.h b/Eigen/src/QR/Tridiagonalization.h
index d580c7e5c..b347b19a6 100755
--- a/Eigen/src/QR/Tridiagonalization.h
+++ b/Eigen/src/QR/Tridiagonalization.h
@@ -163,7 +163,7 @@ void Tridiagonalization<MatrixType>::_compute(MatrixType& matA, CoeffVectorType&
RealScalar beta = ei_sqrt(ei_abs2(v0)+v1norm2);
if (ei_real(v0)>=0.)
beta = -beta;
- matA.col(i).end(n-(i+2)) *= (1./(v0-beta));
+ matA.col(i).end(n-(i+2)) *= (Scalar(1)/(v0-beta));
matA.col(i).coeffRef(i+1) = beta;
Scalar h = (beta - v0) / beta;
// end of the householder transformation
@@ -177,7 +177,7 @@ void Tridiagonalization<MatrixType>::_compute(MatrixType& matA, CoeffVectorType&
* matA.col(i).end(n-i-1));
- hCoeffs.end(n-i-1) += (h * (-0.5) * matA.col(i).end(n-i-1).dot(hCoeffs.end(n-i-1)))
+ hCoeffs.end(n-i-1) += (h * Scalar(-0.5) * matA.col(i).end(n-i-1).dot(hCoeffs.end(n-i-1)))
* matA.col(i).end(n-i-1);
matA.corner(BottomRight,n-i-1,n-i-1).template part<Lower>() =
diff --git a/test/sizeof.cpp b/test/sizeof.cpp
index 4aaeb2df0..6c0485ae1 100644
--- a/test/sizeof.cpp
+++ b/test/sizeof.cpp
@@ -24,7 +24,7 @@
#include "main.h"
-template<typename MatrixType> void verifySizeOf(const MatrixType& m)
+template<typename MatrixType> void verifySizeOf(const MatrixType&)
{
typedef typename MatrixType::Scalar Scalar;
if (MatrixType::RowsAtCompileTime!=Dynamic && MatrixType::ColsAtCompileTime!=Dynamic)