mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-03-19 18:40:38 +08:00
- move: DerivedTraits becomes MatrixBase::Traits
- the static constants are private again in the Derived classes - more documentation and code snippets - new isDiagonal() method
This commit is contained in:
parent
aaf889e72b
commit
84934ea217
@ -66,10 +66,10 @@ template<typename MatrixType, int BlockRows, int BlockCols> class Block
|
||||
|
||||
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Block)
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = BlockRows,
|
||||
ColsAtCompileTime = BlockCols;
|
||||
|
||||
private:
|
||||
const Block& _ref() const { return *this; }
|
||||
int _rows() const { return BlockRows; }
|
||||
int _cols() const { return BlockCols; }
|
||||
|
@ -57,8 +57,8 @@ template<typename NewScalar, typename MatrixType> class Cast : NoOperatorEquals,
|
||||
Cast(const MatRef& matrix) : m_matrix(matrix) {}
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
const Cast& _ref() const { return *this; }
|
||||
int _rows() const { return m_matrix.rows(); }
|
||||
int _cols() const { return m_matrix.cols(); }
|
||||
|
@ -62,10 +62,9 @@ template<typename MatrixType> class Column
|
||||
|
||||
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Column)
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = 1;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = 1;
|
||||
const Column& _ref() const { return *this; }
|
||||
int _rows() const { return m_matrix.rows(); }
|
||||
int _cols() const { return 1; }
|
||||
|
@ -48,10 +48,10 @@ template<typename MatrixType> class Conjugate : NoOperatorEquals,
|
||||
|
||||
Conjugate(const MatRef& matrix) : m_matrix(matrix) {}
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
const Conjugate& _ref() const { return *this; }
|
||||
int _rows() const { return m_matrix.rows(); }
|
||||
int _cols() const { return m_matrix.cols(); }
|
||||
|
@ -50,10 +50,10 @@ template<typename MatrixType> class DiagonalCoeffs
|
||||
|
||||
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(DiagonalCoeffs)
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = 1;
|
||||
|
||||
private:
|
||||
const DiagonalCoeffs& _ref() const { return *this; }
|
||||
int _rows() const { return std::min(m_matrix.rows(), m_matrix.cols()); }
|
||||
int _cols() const { return 1; }
|
||||
|
@ -55,10 +55,10 @@ class DiagonalMatrix : NoOperatorEquals,
|
||||
&& coeffs.size() > 0);
|
||||
}
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = CoeffsVectorType::Traits::SizeAtCompileTime,
|
||||
ColsAtCompileTime = CoeffsVectorType::Traits::SizeAtCompileTime;
|
||||
|
||||
private:
|
||||
|
||||
const DiagonalMatrix& _ref() const { return *this; }
|
||||
int _rows() const { return m_coeffs.size(); }
|
||||
int _cols() const { return m_coeffs.size(); }
|
||||
@ -79,7 +79,7 @@ class DiagonalMatrix : NoOperatorEquals,
|
||||
* Example: \include MatrixBase_asDiagonal.cpp
|
||||
* Output: \verbinclude MatrixBase_asDiagonal.out
|
||||
*
|
||||
* \sa class DiagonalMatrix
|
||||
* \sa class DiagonalMatrix, isDiagonal()
|
||||
**/
|
||||
template<typename Scalar, typename Derived>
|
||||
const DiagonalMatrix<Derived>
|
||||
@ -88,4 +88,32 @@ MatrixBase<Scalar, Derived>::asDiagonal() const
|
||||
return DiagonalMatrix<Derived>(ref());
|
||||
}
|
||||
|
||||
/** \returns true if *this is approximately equal to a diagonal matrix,
|
||||
* within the precision given by \a prec.
|
||||
*
|
||||
* Example: \include MatrixBase_isDiagonal.cpp
|
||||
* Output: \verbinclude MatrixBase_isDiagonal.out
|
||||
*
|
||||
* \sa asDiagonal()
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
bool MatrixBase<Scalar, Derived>::isDiagonal
|
||||
(typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
|
||||
{
|
||||
if(cols() != rows()) return false;
|
||||
RealScalar maxAbsOnDiagonal = static_cast<RealScalar>(-1);
|
||||
for(int j = 0; j < cols(); j++)
|
||||
{
|
||||
RealScalar absOnDiagonal = abs(coeff(j,j));
|
||||
if(absOnDiagonal > maxAbsOnDiagonal) maxAbsOnDiagonal = absOnDiagonal;
|
||||
}
|
||||
for(int j = 0; j < cols(); j++)
|
||||
for(int i = 0; i < j; i++)
|
||||
{
|
||||
if(!Eigen::isMuchSmallerThan(coeff(i, j), maxAbsOnDiagonal, prec)) return false;
|
||||
if(!Eigen::isMuchSmallerThan(coeff(j, i), maxAbsOnDiagonal, prec)) return false;
|
||||
}
|
||||
return true;
|
||||
}
|
||||
|
||||
#endif // EIGEN_DIAGONALMATRIX_H
|
||||
|
@ -41,10 +41,10 @@ template<typename Lhs, typename Rhs> class Difference : NoOperatorEquals,
|
||||
assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
|
||||
}
|
||||
|
||||
static const int RowsAtCompileTime = Lhs::RowsAtCompileTime,
|
||||
ColsAtCompileTime = Rhs::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = Lhs::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = Rhs::Traits::ColsAtCompileTime;
|
||||
|
||||
const Difference& _ref() const { return *this; }
|
||||
int _rows() const { return m_lhs.rows(); }
|
||||
int _cols() const { return m_lhs.cols(); }
|
||||
|
@ -72,10 +72,15 @@ template<typename Scalar, typename Derived>
|
||||
template<typename OtherDerived>
|
||||
Scalar MatrixBase<Scalar, Derived>::dot(const OtherDerived& other) const
|
||||
{
|
||||
assert(Traits::IsVectorAtCompileTime && OtherDerived::Traits::IsVectorAtCompileTime && size() == other.size());
|
||||
assert(Traits::IsVectorAtCompileTime
|
||||
&& OtherDerived::Traits::IsVectorAtCompileTime
|
||||
&& size() == other.size());
|
||||
Scalar res;
|
||||
if(EIGEN_UNROLLED_LOOPS && Traits::SizeAtCompileTime != Dynamic && Traits::SizeAtCompileTime <= 16)
|
||||
DotUnroller<Traits::SizeAtCompileTime-1, Traits::SizeAtCompileTime, Derived, OtherDerived>
|
||||
if(EIGEN_UNROLLED_LOOPS
|
||||
&& Traits::SizeAtCompileTime != Dynamic
|
||||
&& Traits::SizeAtCompileTime <= 16)
|
||||
DotUnroller<Traits::SizeAtCompileTime-1, Traits::SizeAtCompileTime,
|
||||
Derived, OtherDerived>
|
||||
::run(*static_cast<const Derived*>(this), other, res);
|
||||
else
|
||||
{
|
||||
@ -123,25 +128,42 @@ MatrixBase<Scalar, Derived>::normalized() const
|
||||
return (*this) / norm();
|
||||
}
|
||||
|
||||
/** \returns true if *this is approximately orthogonal to \a other,
|
||||
* within the precision given by \a prec.
|
||||
*
|
||||
* Example: \include MatrixBase_isOrtho_vector.cpp
|
||||
* Output: \verbinclude MatrixBase_isOrtho_vector.out
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
template<typename OtherDerived>
|
||||
bool MatrixBase<Scalar, Derived>::isOrtho
|
||||
(const OtherDerived& other,
|
||||
const typename NumTraits<Scalar>::Real& prec = precision<Scalar>()) const
|
||||
typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
|
||||
{
|
||||
return abs2(dot(other)) <= prec * prec * norm2() * other.norm2();
|
||||
}
|
||||
|
||||
/** \returns true if *this is approximately an unitary matrix,
|
||||
* within the precision given by \a prec. In the case where the \a Scalar
|
||||
* type is real numbers, a unitary matrix is an orthogonal matrix, whence the name.
|
||||
*
|
||||
* \note This can be used to check whether a family of vectors forms an orthonormal basis.
|
||||
* Indeed, \c m.isOrtho() returns true if and only if the columns of m form an
|
||||
* orthonormal basis.
|
||||
*
|
||||
* Example: \include MatrixBase_isOrtho_matrix.cpp
|
||||
* Output: \verbinclude MatrixBase_isOrtho_matrix.out
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
bool MatrixBase<Scalar, Derived>::isOrtho
|
||||
(const typename NumTraits<Scalar>::Real& prec = precision<Scalar>()) const
|
||||
(typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
|
||||
{
|
||||
for(int i = 0; i < cols(); i++)
|
||||
{
|
||||
if(!isApprox(col(i).norm2(), static_cast<Scalar>(1)))
|
||||
if(!Eigen::isApprox(col(i).norm2(), static_cast<Scalar>(1), prec))
|
||||
return false;
|
||||
for(int j = 0; j < i; j++)
|
||||
if(!isMuchSmallerThan(col(i).dot(col(j)), static_cast<Scalar>(1)))
|
||||
if(!Eigen::isMuchSmallerThan(col(i).dot(col(j)), static_cast<Scalar>(1), prec))
|
||||
return false;
|
||||
}
|
||||
return true;
|
||||
|
@ -66,11 +66,11 @@ template<typename MatrixType> class DynBlock
|
||||
|
||||
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(DynBlock)
|
||||
|
||||
static const int
|
||||
RowsAtCompileTime = MatrixType::RowsAtCompileTime == 1 ? 1 : Dynamic,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime == 1 ? 1 : Dynamic;
|
||||
|
||||
private:
|
||||
static const int
|
||||
RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime == 1 ? 1 : Dynamic,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime == 1 ? 1 : Dynamic;
|
||||
|
||||
const DynBlock& _ref() const { return *this; }
|
||||
int _rows() const { return m_blockRows; }
|
||||
int _cols() const { return m_blockCols; }
|
||||
|
@ -30,7 +30,7 @@ template<typename Scalar, typename Derived>
|
||||
template<typename OtherDerived>
|
||||
bool MatrixBase<Scalar, Derived>::isApprox(
|
||||
const OtherDerived& other,
|
||||
const typename NumTraits<Scalar>::Real& prec = precision<Scalar>()
|
||||
typename NumTraits<Scalar>::Real prec = precision<Scalar>()
|
||||
) const
|
||||
{
|
||||
assert(rows() == other.rows() && cols() == other.cols());
|
||||
@ -51,7 +51,7 @@ bool MatrixBase<Scalar, Derived>::isApprox(
|
||||
template<typename Scalar, typename Derived>
|
||||
bool MatrixBase<Scalar, Derived>::isMuchSmallerThan(
|
||||
const typename NumTraits<Scalar>::Real& other,
|
||||
const typename NumTraits<Scalar>::Real& prec = precision<Scalar>()
|
||||
typename NumTraits<Scalar>::Real prec = precision<Scalar>()
|
||||
) const
|
||||
{
|
||||
if(Traits::IsVectorAtCompileTime)
|
||||
@ -71,7 +71,7 @@ template<typename Scalar, typename Derived>
|
||||
template<typename OtherDerived>
|
||||
bool MatrixBase<Scalar, Derived>::isMuchSmallerThan(
|
||||
const MatrixBase<Scalar, OtherDerived>& other,
|
||||
const typename NumTraits<Scalar>::Real& prec = precision<Scalar>()
|
||||
typename NumTraits<Scalar>::Real prec = precision<Scalar>()
|
||||
) const
|
||||
{
|
||||
assert(rows() == other.rows() && cols() == other.cols());
|
||||
|
@ -26,6 +26,12 @@
|
||||
#ifndef EIGEN_IDENTITY_H
|
||||
#define EIGEN_IDENTITY_H
|
||||
|
||||
/** \class Identity
|
||||
*
|
||||
* \brief Expression of the identity matrix of some size.
|
||||
*
|
||||
* \sa MatrixBase::identity(int)
|
||||
*/
|
||||
template<typename MatrixType> class Identity : NoOperatorEquals,
|
||||
public MatrixBase<typename MatrixType::Scalar, Identity<MatrixType> >
|
||||
{
|
||||
@ -33,15 +39,15 @@ template<typename MatrixType> class Identity : NoOperatorEquals,
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
friend class MatrixBase<Scalar, Identity<MatrixType> >;
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
Identity(int rows) : m_rows(rows)
|
||||
{
|
||||
assert(rows > 0 && RowsAtCompileTime == ColsAtCompileTime);
|
||||
}
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
const Identity& _ref() const { return *this; }
|
||||
int _rows() const { return m_rows; }
|
||||
int _cols() const { return m_rows; }
|
||||
@ -55,22 +61,46 @@ template<typename MatrixType> class Identity : NoOperatorEquals,
|
||||
int m_rows;
|
||||
};
|
||||
|
||||
/** \returns an expression of the identity matrix of given type and size.
|
||||
*
|
||||
* \param rows The number of rows of the identity matrix to return. If *this has
|
||||
* fixed size, that size is used as the default argument for \a rows
|
||||
* and is then the only allowed value. If *this has dynamic size,
|
||||
* you can use any positive value for \a rows.
|
||||
*
|
||||
* \note An identity matrix is a square matrix, so it is required that the type of *this
|
||||
* allows being a square matrix.
|
||||
*
|
||||
* Example: \include MatrixBase_identity_int.cpp
|
||||
* Output: \verbinclude MatrixBase_identity_int.out
|
||||
*
|
||||
* \sa class Identity, isIdentity()
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
const Identity<Derived> MatrixBase<Scalar, Derived>::identity(int rows)
|
||||
{
|
||||
return Identity<Derived>(rows);
|
||||
}
|
||||
|
||||
/** \returns true if *this is approximately equal to the identity matrix,
|
||||
* within the precision given by \a prec.
|
||||
*
|
||||
* Example: \include MatrixBase_isIdentity.cpp
|
||||
* Output: \verbinclude MatrixBase_isIdentity.out
|
||||
*
|
||||
* \sa class Identity, identity(int)
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
bool MatrixBase<Scalar, Derived>::isIdentity
|
||||
(const typename NumTraits<Scalar>::Real& prec = precision<Scalar>()) const
|
||||
(typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
|
||||
{
|
||||
for(int j = 0; j < col(); j++)
|
||||
if(cols() != rows()) return false;
|
||||
for(int j = 0; j < cols(); j++)
|
||||
{
|
||||
if(!isApprox(coeff(j, j), static_cast<Scalar>(1)))
|
||||
if(!Eigen::isApprox(coeff(j, j), static_cast<Scalar>(1), prec))
|
||||
return false;
|
||||
for(int i = 0; i < j; i++)
|
||||
if(!isMuchSmallerThan(coeff(i, j), static_cast<Scalar>(1)))
|
||||
if(!Eigen::isMuchSmallerThan(coeff(i, j), static_cast<Scalar>(1), prec))
|
||||
return false;
|
||||
}
|
||||
return true;
|
||||
|
@ -46,18 +46,18 @@ template<typename MatrixType> class Map
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
friend class MatrixBase<Scalar, Map<MatrixType> >;
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
static const MatrixStorageOrder Order = MatrixType::Order;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
static const MatrixStorageOrder _Order = MatrixType::StorageOrder;
|
||||
|
||||
const Map& _ref() const { return *this; }
|
||||
int _rows() const { return m_rows; }
|
||||
int _cols() const { return m_cols; }
|
||||
|
||||
const Scalar& _coeff(int row, int col) const
|
||||
{
|
||||
if(Order == ColumnMajor)
|
||||
if(_Order == ColumnMajor)
|
||||
return m_data[row + col * m_rows];
|
||||
else // RowMajor
|
||||
return m_data[col + row * m_cols];
|
||||
@ -65,7 +65,7 @@ template<typename MatrixType> class Map
|
||||
|
||||
Scalar& _coeffRef(int row, int col)
|
||||
{
|
||||
if(Order == ColumnMajor)
|
||||
if(_Order == ColumnMajor)
|
||||
return const_cast<Scalar*>(m_data)[row + col * m_rows];
|
||||
else // RowMajor
|
||||
return const_cast<Scalar*>(m_data)[col + row * m_cols];
|
||||
|
@ -144,6 +144,10 @@ inline bool isMuchSmallerThan(const std::complex<float>& a, const std::complex<f
|
||||
{
|
||||
return abs2(a) <= abs2(b) * prec * prec;
|
||||
}
|
||||
inline bool isMuchSmallerThan(const std::complex<float>& a, float b, float prec = precision<float>())
|
||||
{
|
||||
return abs2(a) <= abs2(b) * prec * prec;
|
||||
}
|
||||
inline bool isApprox(const std::complex<float>& a, const std::complex<float>& b, float prec = precision<float>())
|
||||
{
|
||||
return isApprox(std::real(a), std::real(b), prec)
|
||||
@ -165,6 +169,10 @@ inline bool isMuchSmallerThan(const std::complex<double>& a, const std::complex<
|
||||
{
|
||||
return abs2(a) <= abs2(b) * prec * prec;
|
||||
}
|
||||
inline bool isMuchSmallerThan(const std::complex<double>& a, double b, double prec = precision<double>())
|
||||
{
|
||||
return abs2(a) <= abs2(b) * prec * prec;
|
||||
}
|
||||
inline bool isApprox(const std::complex<double>& a, const std::complex<double>& b, double prec = precision<double>())
|
||||
{
|
||||
return isApprox(std::real(a), std::real(b), prec)
|
||||
|
@ -79,6 +79,8 @@ class Matrix : public MatrixBase<_Scalar, Matrix<_Scalar, _Rows, _Cols, _Storage
|
||||
{
|
||||
public:
|
||||
friend class MatrixBase<_Scalar, Matrix>;
|
||||
friend class Map<Matrix>;
|
||||
|
||||
typedef MatrixBase<_Scalar, Matrix> Base;
|
||||
typedef MatrixStorage<_Scalar, _Rows, _Cols> Storage;
|
||||
typedef _Scalar Scalar;
|
||||
@ -91,16 +93,15 @@ class Matrix : public MatrixBase<_Scalar, Matrix<_Scalar, _Rows, _Cols, _Storage
|
||||
Scalar* data()
|
||||
{ return Storage::m_data; }
|
||||
|
||||
static const MatrixStorageOrder Order = _StorageOrder;
|
||||
static const int RowsAtCompileTime = _Rows, ColsAtCompileTime = _Cols;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = _Rows, ColsAtCompileTime = _Cols;
|
||||
static const MatrixStorageOrder StorageOrder = _StorageOrder;
|
||||
|
||||
Ref _ref() const { return Ref(*this); }
|
||||
|
||||
const Scalar& _coeff(int row, int col) const
|
||||
{
|
||||
if(Order == ColumnMajor)
|
||||
if(_StorageOrder == ColumnMajor)
|
||||
return (Storage::m_data)[row + col * Storage::_rows()];
|
||||
else // RowMajor
|
||||
return (Storage::m_data)[col + row * Storage::_cols()];
|
||||
@ -108,7 +109,7 @@ class Matrix : public MatrixBase<_Scalar, Matrix<_Scalar, _Rows, _Cols, _Storage
|
||||
|
||||
Scalar& _coeffRef(int row, int col)
|
||||
{
|
||||
if(Order == ColumnMajor)
|
||||
if(_StorageOrder == ColumnMajor)
|
||||
return (Storage::m_data)[row + col * Storage::_rows()];
|
||||
else // RowMajor
|
||||
return (Storage::m_data)[col + row * Storage::_cols()];
|
||||
|
@ -26,37 +26,6 @@
|
||||
#ifndef EIGEN_MATRIXBASE_H
|
||||
#define EIGEN_MATRIXBASE_H
|
||||
|
||||
#include "Util.h"
|
||||
|
||||
template <typename Derived>
|
||||
struct DerivedTraits
|
||||
{
|
||||
/** The number of rows at compile-time. This is just a copy of the value provided
|
||||
* by the \a Derived type. If a value is not known at compile-time,
|
||||
* it is set to the \a Dynamic constant.
|
||||
* \sa rows(), cols(), ColsAtCompileTime, SizeAtCompileTime */
|
||||
static const int RowsAtCompileTime = Derived::RowsAtCompileTime;
|
||||
|
||||
/** The number of columns at compile-time. This is just a copy of the value provided
|
||||
* by the \a Derived type. If a value is not known at compile-time,
|
||||
* it is set to the \a Dynamic constant.
|
||||
* \sa rows(), cols(), RowsAtCompileTime, SizeAtCompileTime */
|
||||
static const int ColsAtCompileTime = Derived::ColsAtCompileTime;
|
||||
|
||||
/** This is equal to the number of coefficients, i.e. the number of
|
||||
* rows times the number of columns, or to \a Dynamic if this is not
|
||||
* known at compile-time. \sa RowsAtCompileTime, ColsAtCompileTime */
|
||||
static const int SizeAtCompileTime
|
||||
= Derived::RowsAtCompileTime == Dynamic || Derived::ColsAtCompileTime == Dynamic
|
||||
? Dynamic : Derived::RowsAtCompileTime * Derived::ColsAtCompileTime;
|
||||
|
||||
/** This is set to true if either the number of rows or the number of
|
||||
* columns is known at compile-time to be equal to 1. Indeed, in that case,
|
||||
* we are dealing with a column-vector (if there is only one column) or with
|
||||
* a row-vector (if there is only one row). */
|
||||
static const bool IsVectorAtCompileTime = Derived::RowsAtCompileTime == 1 || Derived::ColsAtCompileTime == 1;
|
||||
};
|
||||
|
||||
/** \class MatrixBase
|
||||
*
|
||||
* \brief Base class for all matrices, vectors, and expressions
|
||||
@ -88,7 +57,37 @@ struct DerivedTraits
|
||||
template<typename Scalar, typename Derived> class MatrixBase
|
||||
{
|
||||
public:
|
||||
typedef DerivedTraits<Derived> Traits;
|
||||
|
||||
/** \brief Some traits provided by the Derived type.
|
||||
* Grouping these in a nested subclass is what was needed for ICC compatibility. */
|
||||
struct Traits
|
||||
{
|
||||
/** The number of rows at compile-time. This is just a copy of the value provided
|
||||
* by the \a Derived type. If a value is not known at compile-time,
|
||||
* it is set to the \a Dynamic constant.
|
||||
* \sa MatrixBase::rows(), MatrixBase::cols(), ColsAtCompileTime, SizeAtCompileTime */
|
||||
static const int RowsAtCompileTime = Derived::RowsAtCompileTime;
|
||||
|
||||
/** The number of columns at compile-time. This is just a copy of the value provided
|
||||
* by the \a Derived type. If a value is not known at compile-time,
|
||||
* it is set to the \a Dynamic constant.
|
||||
* \sa MatrixBase::rows(), MatrixBase::cols(), RowsAtCompileTime, SizeAtCompileTime */
|
||||
static const int ColsAtCompileTime = Derived::ColsAtCompileTime;
|
||||
|
||||
/** This is equal to the number of coefficients, i.e. the number of
|
||||
* rows times the number of columns, or to \a Dynamic if this is not
|
||||
* known at compile-time. \sa RowsAtCompileTime, ColsAtCompileTime */
|
||||
static const int SizeAtCompileTime
|
||||
= Derived::RowsAtCompileTime == Dynamic || Derived::ColsAtCompileTime == Dynamic
|
||||
? Dynamic : Derived::RowsAtCompileTime * Derived::ColsAtCompileTime;
|
||||
|
||||
/** This is set to true if either the number of rows or the number of
|
||||
* columns is known at compile-time to be equal to 1. Indeed, in that case,
|
||||
* we are dealing with a column-vector (if there is only one column) or with
|
||||
* a row-vector (if there is only one row). */
|
||||
static const bool IsVectorAtCompileTime
|
||||
= Derived::RowsAtCompileTime == 1 || Derived::ColsAtCompileTime == 1;
|
||||
};
|
||||
|
||||
/** This is the "reference type" used to pass objects of type MatrixBase as arguments
|
||||
* to functions. If this MatrixBase type represents an expression, then \a Ref
|
||||
@ -104,17 +103,17 @@ template<typename Scalar, typename Derived> class MatrixBase
|
||||
* \a Scalar is \a std::complex<T> then RealScalar is \a T. */
|
||||
typedef typename NumTraits<Scalar>::Real RealScalar;
|
||||
|
||||
/** \returns the number of rows. \sa cols(), RowsAtCompileTime */
|
||||
/** \returns the number of rows. \sa cols(), DerivedTraits::RowsAtCompileTime */
|
||||
int rows() const { return static_cast<const Derived *>(this)->_rows(); }
|
||||
/** \returns the number of columns. \sa row(), ColsAtCompileTime*/
|
||||
/** \returns the number of columns. \sa row(), DerivedTraits::ColsAtCompileTime*/
|
||||
int cols() const { return static_cast<const Derived *>(this)->_cols(); }
|
||||
/** \returns the number of coefficients, which is \a rows()*cols().
|
||||
* \sa rows(), cols(), SizeAtCompileTime. */
|
||||
* \sa rows(), cols(), DerivedTraits::SizeAtCompileTime. */
|
||||
int size() const { return rows() * cols(); }
|
||||
/** \returns true if either the number of rows or the number of columns is equal to 1.
|
||||
* In other words, this function returns
|
||||
* \code rows()==1 || cols()==1 \endcode
|
||||
* \sa rows(), cols(), IsVectorAtCompileTime. */
|
||||
* \sa rows(), cols(), DerivedTraits::IsVectorAtCompileTime. */
|
||||
bool isVector() const { return rows()==1 || cols()==1; }
|
||||
/** \returns a Ref to *this. \sa Ref */
|
||||
Ref ref() const
|
||||
@ -164,9 +163,8 @@ template<typename Scalar, typename Derived> class MatrixBase
|
||||
RealScalar norm() const;
|
||||
const ScalarMultiple<RealScalar, Derived> normalized() const;
|
||||
template<typename OtherDerived>
|
||||
bool isOrtho(const OtherDerived& other,
|
||||
const typename NumTraits<Scalar>::Real& prec) const;
|
||||
bool isOrtho(const typename NumTraits<Scalar>::Real& prec) const;
|
||||
bool isOrtho(const OtherDerived& other, RealScalar prec) const;
|
||||
bool isOrtho(RealScalar prec) const;
|
||||
|
||||
static const Eval<Random<Derived> > random(int rows, int cols);
|
||||
static const Eval<Random<Derived> > random(int size);
|
||||
@ -179,9 +177,10 @@ template<typename Scalar, typename Derived> class MatrixBase
|
||||
static const Ones<Derived> ones();
|
||||
static const Identity<Derived> identity(int rows = Derived::RowsAtCompileTime);
|
||||
|
||||
bool isZero(const typename NumTraits<Scalar>::Real& prec) const;
|
||||
bool isOnes(const typename NumTraits<Scalar>::Real& prec) const;
|
||||
bool isIdentity(const typename NumTraits<Scalar>::Real& prec) const;
|
||||
bool isZero(RealScalar prec) const;
|
||||
bool isOnes(RealScalar prec) const;
|
||||
bool isIdentity(RealScalar prec) const;
|
||||
bool isDiagonal(RealScalar prec) const;
|
||||
|
||||
const DiagonalMatrix<Derived> asDiagonal() const;
|
||||
|
||||
@ -189,19 +188,12 @@ template<typename Scalar, typename Derived> class MatrixBase
|
||||
const DiagonalCoeffs<Derived> diagonal() const;
|
||||
|
||||
template<typename OtherDerived>
|
||||
bool isApprox(
|
||||
const OtherDerived& other,
|
||||
const typename NumTraits<Scalar>::Real& prec
|
||||
) const;
|
||||
bool isMuchSmallerThan(
|
||||
const typename NumTraits<Scalar>::Real& other,
|
||||
const typename NumTraits<Scalar>::Real& prec
|
||||
) const;
|
||||
bool isApprox(const OtherDerived& other, RealScalar prec) const;
|
||||
bool isMuchSmallerThan
|
||||
(const typename NumTraits<Scalar>::Real& other, RealScalar prec) const;
|
||||
template<typename OtherDerived>
|
||||
bool isMuchSmallerThan(
|
||||
const MatrixBase<Scalar, OtherDerived>& other,
|
||||
const typename NumTraits<Scalar>::Real& prec
|
||||
) const;
|
||||
bool isMuchSmallerThan
|
||||
(const MatrixBase<Scalar, OtherDerived>& other, RealScalar prec) const;
|
||||
|
||||
template<typename OtherDerived>
|
||||
const Product<Derived, OtherDerived>
|
||||
|
@ -38,10 +38,10 @@ template<typename MatrixType> class MatrixRef
|
||||
|
||||
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(MatrixRef)
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
int _rows() const { return m_matrix.rows(); }
|
||||
int _cols() const { return m_matrix.cols(); }
|
||||
|
||||
|
@ -56,13 +56,13 @@ template<typename MatrixType> class Minor
|
||||
|
||||
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Minor)
|
||||
|
||||
static const int
|
||||
RowsAtCompileTime = (MatrixType::RowsAtCompileTime != Dynamic) ?
|
||||
MatrixType::RowsAtCompileTime - 1 : Dynamic,
|
||||
ColsAtCompileTime = (MatrixType::ColsAtCompileTime != Dynamic) ?
|
||||
MatrixType::ColsAtCompileTime - 1 : Dynamic;
|
||||
|
||||
private:
|
||||
static const int
|
||||
RowsAtCompileTime = (MatrixType::Traits::RowsAtCompileTime != Dynamic) ?
|
||||
MatrixType::Traits::RowsAtCompileTime - 1 : Dynamic,
|
||||
ColsAtCompileTime = (MatrixType::Traits::ColsAtCompileTime != Dynamic) ?
|
||||
MatrixType::Traits::ColsAtCompileTime - 1 : Dynamic;
|
||||
|
||||
const Minor& _ref() const { return *this; }
|
||||
int _rows() const { return m_matrix.rows() - 1; }
|
||||
int _cols() const { return m_matrix.cols() - 1; }
|
||||
|
@ -39,10 +39,9 @@ template<typename MatrixType> class Ones : NoOperatorEquals,
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
friend class MatrixBase<Scalar, Ones<MatrixType> >;
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
const Ones& _ref() const { return *this; }
|
||||
int _rows() const { return m_rows; }
|
||||
@ -78,7 +77,7 @@ template<typename MatrixType> class Ones : NoOperatorEquals,
|
||||
* Example: \include MatrixBase_ones_int_int.cpp
|
||||
* Output: \verbinclude MatrixBase_ones_int_int.out
|
||||
*
|
||||
* \sa ones(), ones(int)
|
||||
* \sa ones(), ones(int), isOnes(), class Ones
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
const Ones<Derived> MatrixBase<Scalar, Derived>::ones(int rows, int cols)
|
||||
@ -100,7 +99,7 @@ const Ones<Derived> MatrixBase<Scalar, Derived>::ones(int rows, int cols)
|
||||
* Example: \include MatrixBase_ones_int.cpp
|
||||
* Output: \verbinclude MatrixBase_ones_int.out
|
||||
*
|
||||
* \sa ones(), ones(int,int)
|
||||
* \sa ones(), ones(int,int), isOnes(), class Ones
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
const Ones<Derived> MatrixBase<Scalar, Derived>::ones(int size)
|
||||
@ -118,7 +117,7 @@ const Ones<Derived> MatrixBase<Scalar, Derived>::ones(int size)
|
||||
* Example: \include MatrixBase_ones.cpp
|
||||
* Output: \verbinclude MatrixBase_ones.out
|
||||
*
|
||||
* \sa ones(int), ones(int,int)
|
||||
* \sa ones(int), ones(int,int), isOnes(), class Ones
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
const Ones<Derived> MatrixBase<Scalar, Derived>::ones()
|
||||
@ -126,13 +125,21 @@ const Ones<Derived> MatrixBase<Scalar, Derived>::ones()
|
||||
return Ones<Derived>(Traits::RowsAtCompileTime, Traits::ColsAtCompileTime);
|
||||
}
|
||||
|
||||
/** \returns true if *this is approximately equal to the matrix where all coefficients
|
||||
* are equal to 1, within the precision given by \a prec.
|
||||
*
|
||||
* Example: \include MatrixBase_isOnes.cpp
|
||||
* Output: \verbinclude MatrixBase_isOnes.out
|
||||
*
|
||||
* \sa class Ones, ones()
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
bool MatrixBase<Scalar, Derived>::isOnes
|
||||
(const typename NumTraits<Scalar>::Real& prec = precision<Scalar>()) const
|
||||
(typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
|
||||
{
|
||||
for(int j = 0; j < col(); j++)
|
||||
for(int i = 0; i < row(); i++)
|
||||
if(!isApprox(coeff(i, j), static_cast<Scalar>(1)))
|
||||
for(int j = 0; j < cols(); j++)
|
||||
for(int i = 0; i < rows(); i++)
|
||||
if(!Eigen::isApprox(coeff(i, j), static_cast<Scalar>(1), prec))
|
||||
return false;
|
||||
return true;
|
||||
}
|
||||
|
@ -30,8 +30,8 @@
|
||||
template<typename Derived1, typename Derived2, int UnrollCount>
|
||||
struct MatrixOperatorEqualsUnroller
|
||||
{
|
||||
static const int col = (UnrollCount-1) / Derived1::RowsAtCompileTime;
|
||||
static const int row = (UnrollCount-1) % Derived1::RowsAtCompileTime;
|
||||
static const int col = (UnrollCount-1) / Derived1::Traits::RowsAtCompileTime;
|
||||
static const int row = (UnrollCount-1) % Derived1::Traits::RowsAtCompileTime;
|
||||
|
||||
static void run(Derived1 &dst, const Derived2 &src)
|
||||
{
|
||||
|
@ -36,10 +36,10 @@ template<typename MatrixType> class Opposite : NoOperatorEquals,
|
||||
|
||||
Opposite(const MatRef& matrix) : m_matrix(matrix) {}
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
const Opposite& _ref() const { return *this; }
|
||||
int _rows() const { return m_matrix.rows(); }
|
||||
int _cols() const { return m_matrix.cols(); }
|
||||
|
@ -75,10 +75,9 @@ template<typename Lhs, typename Rhs> class Product : NoOperatorEquals,
|
||||
assert(lhs.cols() == rhs.rows());
|
||||
}
|
||||
|
||||
static const int RowsAtCompileTime = Lhs::RowsAtCompileTime,
|
||||
ColsAtCompileTime = Rhs::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = Lhs::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = Rhs::Traits::ColsAtCompileTime;
|
||||
|
||||
const Product& _ref() const { return *this; }
|
||||
int _rows() const { return m_lhs.rows(); }
|
||||
@ -88,8 +87,10 @@ template<typename Lhs, typename Rhs> class Product : NoOperatorEquals,
|
||||
{
|
||||
Scalar res;
|
||||
if(EIGEN_UNROLLED_LOOPS
|
||||
&& Lhs::ColsAtCompileTime != Dynamic && Lhs::ColsAtCompileTime <= 16)
|
||||
ProductUnroller<Lhs::ColsAtCompileTime-1, Lhs::ColsAtCompileTime, LhsRef, RhsRef>
|
||||
&& Lhs::Traits::ColsAtCompileTime != Dynamic
|
||||
&& Lhs::Traits::ColsAtCompileTime <= 16)
|
||||
ProductUnroller<Lhs::Traits::ColsAtCompileTime-1,
|
||||
Lhs::Traits::ColsAtCompileTime, LhsRef, RhsRef>
|
||||
::run(row, col, m_lhs, m_rhs, res);
|
||||
else
|
||||
{
|
||||
|
@ -39,10 +39,9 @@ template<typename MatrixType> class Random : NoOperatorEquals,
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
friend class MatrixBase<Scalar, Random<MatrixType> >;
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
const Random& _ref() const { return *this; }
|
||||
int _rows() const { return m_rows; }
|
||||
@ -50,7 +49,7 @@ template<typename MatrixType> class Random : NoOperatorEquals,
|
||||
|
||||
Scalar _coeff(int, int) const
|
||||
{
|
||||
return random<Scalar>();
|
||||
return Eigen::random<Scalar>();
|
||||
}
|
||||
|
||||
public:
|
||||
|
@ -68,10 +68,10 @@ template<typename MatrixType> class Row
|
||||
|
||||
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Row)
|
||||
|
||||
static const int RowsAtCompileTime = 1,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = 1,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
const Row& _ref() const { return *this; }
|
||||
|
||||
int _rows() const { return 1; }
|
||||
|
@ -37,10 +37,10 @@ template<typename FactorType, typename MatrixType> class ScalarMultiple : NoOper
|
||||
ScalarMultiple(const MatRef& matrix, FactorType factor)
|
||||
: m_matrix(matrix), m_factor(factor) {}
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
const ScalarMultiple& _ref() const { return *this; }
|
||||
int _rows() const { return m_matrix.rows(); }
|
||||
int _cols() const { return m_matrix.cols(); }
|
||||
|
@ -41,10 +41,10 @@ template<typename Lhs, typename Rhs> class Sum : NoOperatorEquals,
|
||||
assert(lhs.rows() == rhs.rows() && lhs.cols() == rhs.cols());
|
||||
}
|
||||
|
||||
static const int RowsAtCompileTime = Lhs::RowsAtCompileTime,
|
||||
ColsAtCompileTime = Rhs::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = Lhs::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = Rhs::Traits::ColsAtCompileTime;
|
||||
|
||||
const Sum& _ref() const { return *this; }
|
||||
int _rows() const { return m_lhs.rows(); }
|
||||
int _cols() const { return m_lhs.cols(); }
|
||||
|
@ -50,10 +50,10 @@ template<typename MatrixType> class Transpose
|
||||
|
||||
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Transpose)
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::ColsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::RowsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::ColsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::RowsAtCompileTime;
|
||||
|
||||
const Transpose& _ref() const { return *this; }
|
||||
int _rows() const { return m_matrix.cols(); }
|
||||
int _cols() const { return m_matrix.rows(); }
|
||||
|
@ -39,10 +39,10 @@ template<typename MatrixType> class Zero : NoOperatorEquals,
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
friend class MatrixBase<Scalar, Zero<MatrixType> >;
|
||||
|
||||
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
|
||||
|
||||
private:
|
||||
static const int RowsAtCompileTime = MatrixType::Traits::RowsAtCompileTime,
|
||||
ColsAtCompileTime = MatrixType::Traits::ColsAtCompileTime;
|
||||
|
||||
const Zero& _ref() const { return *this; }
|
||||
int _rows() const { return m_rows; }
|
||||
int _cols() const { return m_cols; }
|
||||
@ -125,13 +125,21 @@ const Zero<Derived> MatrixBase<Scalar, Derived>::zero()
|
||||
return Zero<Derived>(Traits::RowsAtCompileTime, Traits::ColsAtCompileTime);
|
||||
}
|
||||
|
||||
/** \returns true if *this is approximately equal to the zero matrix,
|
||||
* within the precision given by \a prec.
|
||||
*
|
||||
* Example: \include MatrixBase_isZero.cpp
|
||||
* Output: \verbinclude MatrixBase_isZero.out
|
||||
*
|
||||
* \sa class Zero, zero()
|
||||
*/
|
||||
template<typename Scalar, typename Derived>
|
||||
bool MatrixBase<Scalar, Derived>::isZero
|
||||
(const typename NumTraits<Scalar>::Real& prec = precision<Scalar>()) const
|
||||
(typename NumTraits<Scalar>::Real prec = precision<Scalar>()) const
|
||||
{
|
||||
for(int j = 0; j < col(); j++)
|
||||
for(int i = 0; i < row(); i++)
|
||||
if(!isMuchSmallerThan(coeff(i, j), static_cast<Scalar>(1)))
|
||||
for(int j = 0; j < cols(); j++)
|
||||
for(int i = 0; i < rows(); i++)
|
||||
if(!Eigen::isMuchSmallerThan(coeff(i, j), static_cast<Scalar>(1), prec))
|
||||
return false;
|
||||
return true;
|
||||
}
|
||||
|
@ -37,7 +37,7 @@ SEPARATE_MEMBER_PAGES = NO
|
||||
TAB_SIZE = 8
|
||||
ALIASES = \
|
||||
"only_for_vectors=This is only for vectors (either row-vectors or column-vectors), \
|
||||
as determined by \link MatrixBase::IsVectorAtCompileTime \endlink."
|
||||
i.e. matrices which are known at compile-time to have either one row or one column."
|
||||
OPTIMIZE_OUTPUT_FOR_C = NO
|
||||
OPTIMIZE_OUTPUT_JAVA = NO
|
||||
BUILTIN_STL_SUPPORT = NO
|
||||
|
6
doc/snippets/MatrixBase_isDiagonal.cpp
Normal file
6
doc/snippets/MatrixBase_isDiagonal.cpp
Normal file
@ -0,0 +1,6 @@
|
||||
Matrix3d m = 10000 * Matrix3d::identity();
|
||||
m(0,2) = 1;
|
||||
cout << "Here's the matrix m:" << endl << m << endl;
|
||||
cout << "m.isDiagonal() returns: " << m.isDiagonal() << endl;
|
||||
cout << "m.isDiagonal(1e-3) returns: " << m.isDiagonal(1e-3) << endl;
|
||||
|
5
doc/snippets/MatrixBase_isIdentity.cpp
Normal file
5
doc/snippets/MatrixBase_isIdentity.cpp
Normal file
@ -0,0 +1,5 @@
|
||||
Matrix3d m = Matrix3d::identity();
|
||||
m(0,2) = 1e-4;
|
||||
cout << "Here's the matrix m:" << endl << m << endl;
|
||||
cout << "m.isIdentity() returns: " << m.isIdentity() << endl;
|
||||
cout << "m.isIdentity(1e-3) returns: " << m.isIdentity(1e-3) << endl;
|
5
doc/snippets/MatrixBase_isOnes.cpp
Normal file
5
doc/snippets/MatrixBase_isOnes.cpp
Normal file
@ -0,0 +1,5 @@
|
||||
Matrix3d m = Matrix3d::ones();
|
||||
m(0,2) += 1e-4;
|
||||
cout << "Here's the matrix m:" << endl << m << endl;
|
||||
cout << "m.isOnes() returns: " << m.isOnes() << endl;
|
||||
cout << "m.isOnes(1e-3) returns: " << m.isOnes(1e-3) << endl;
|
5
doc/snippets/MatrixBase_isOrtho_matrix.cpp
Normal file
5
doc/snippets/MatrixBase_isOrtho_matrix.cpp
Normal file
@ -0,0 +1,5 @@
|
||||
Matrix3d m = Matrix3d::identity();
|
||||
m(0,2) = 1e-4;
|
||||
cout << "Here's the matrix m:" << endl << m << endl;
|
||||
cout << "m.isOrtho() returns: " << m.isOrtho() << endl;
|
||||
cout << "m.isOrtho(1e-3) returns: " << m.isOrtho(1e-3) << endl;
|
6
doc/snippets/MatrixBase_isOrtho_vector.cpp
Normal file
6
doc/snippets/MatrixBase_isOrtho_vector.cpp
Normal file
@ -0,0 +1,6 @@
|
||||
Vector3d v(1,0,0);
|
||||
Vector3d w(1e-4,0,1);
|
||||
cout << "Here's the vector v:" << endl << v << endl;
|
||||
cout << "Here's the vector w:" << endl << w << endl;
|
||||
cout << "v.isOrtho(w) returns: " << v.isOrtho(w) << endl;
|
||||
cout << "v.isOrtho(w,1e-3) returns: " << v.isOrtho(w,1e-3) << endl;
|
5
doc/snippets/MatrixBase_isZero.cpp
Normal file
5
doc/snippets/MatrixBase_isZero.cpp
Normal file
@ -0,0 +1,5 @@
|
||||
Matrix3d m = Matrix3d::zero();
|
||||
m(0,2) = 1e-4;
|
||||
cout << "Here's the matrix m:" << endl << m << endl;
|
||||
cout << "m.isZero() returns: " << m.isZero() << endl;
|
||||
cout << "m.isZero(1e-3) returns: " << m.isZero(1e-3) << endl;
|
@ -34,7 +34,7 @@ template<typename MatrixType> void adjoint(const MatrixType& m)
|
||||
*/
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, 1> VectorType;
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
@ -42,9 +42,9 @@ template<typename MatrixType> void adjoint(const MatrixType& m)
|
||||
m2 = MatrixType::random(rows, cols),
|
||||
m3(rows, cols),
|
||||
mzero = MatrixType::zero(rows, cols),
|
||||
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
identity = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::identity(rows),
|
||||
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
square = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::random(rows, rows);
|
||||
VectorType v1 = VectorType::random(rows),
|
||||
v2 = VectorType::random(rows),
|
||||
|
@ -30,7 +30,7 @@ namespace Eigen {
|
||||
template<typename MatrixType> void basicStuff(const MatrixType& m)
|
||||
{
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, 1> VectorType;
|
||||
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
@ -41,9 +41,9 @@ template<typename MatrixType> void basicStuff(const MatrixType& m)
|
||||
m2 = MatrixType::random(rows, cols),
|
||||
m3(rows, cols),
|
||||
mzero = MatrixType::zero(rows, cols),
|
||||
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
identity = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::identity(rows),
|
||||
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
square = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::random(rows, rows);
|
||||
VectorType v1 = VectorType::random(rows),
|
||||
v2 = VectorType::random(rows),
|
||||
@ -75,8 +75,8 @@ template<typename MatrixType> void basicStuff(const MatrixType& m)
|
||||
|
||||
// now test copying a row-vector into a (column-)vector and conversely.
|
||||
square.col(r) = square.row(r).eval();
|
||||
Matrix<Scalar, 1, MatrixType::RowsAtCompileTime> rv(rows);
|
||||
Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> cv(rows);
|
||||
Matrix<Scalar, 1, MatrixType::Traits::RowsAtCompileTime> rv(rows);
|
||||
Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, 1> cv(rows);
|
||||
rv = square.col(r);
|
||||
cv = square.row(r);
|
||||
VERIFY_IS_APPROX(rv, cv.transpose());
|
||||
|
@ -34,7 +34,7 @@ template<typename MatrixType> void linearStructure(const MatrixType& m)
|
||||
*/
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, 1> VectorType;
|
||||
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
@ -45,9 +45,9 @@ template<typename MatrixType> void linearStructure(const MatrixType& m)
|
||||
m2 = MatrixType::random(rows, cols),
|
||||
m3(rows, cols),
|
||||
mzero = MatrixType::zero(rows, cols),
|
||||
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
identity = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::identity(rows),
|
||||
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
square = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::random(rows, rows);
|
||||
VectorType v1 = VectorType::random(rows),
|
||||
v2 = VectorType::random(rows),
|
||||
|
@ -34,8 +34,8 @@ template<typename MatrixType> void miscMatrices(const MatrixType& m)
|
||||
*/
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
|
||||
typedef Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, 1, MatrixType::Traits::ColsAtCompileTime> RowVectorType;
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
@ -45,7 +45,7 @@ template<typename MatrixType> void miscMatrices(const MatrixType& m)
|
||||
VERIFY_IS_APPROX(m1(r,c), static_cast<Scalar>(1));
|
||||
VectorType v1 = VectorType::random(rows);
|
||||
v1[0];
|
||||
Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
square = v1.asDiagonal();
|
||||
if(r==r2) VERIFY_IS_APPROX(square(r,r2), v1[r]);
|
||||
else VERIFY_IS_MUCH_SMALLER_THAN(square(r,r2), static_cast<Scalar>(1));
|
||||
|
@ -34,7 +34,7 @@ template<typename MatrixType> void product(const MatrixType& m)
|
||||
*/
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, 1> VectorType;
|
||||
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
@ -45,9 +45,9 @@ template<typename MatrixType> void product(const MatrixType& m)
|
||||
m2 = MatrixType::random(rows, cols),
|
||||
m3(rows, cols),
|
||||
mzero = MatrixType::zero(rows, cols),
|
||||
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
identity = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::identity(rows),
|
||||
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
square = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::random(rows, rows);
|
||||
VectorType v1 = VectorType::random(rows),
|
||||
v2 = VectorType::random(rows),
|
||||
|
@ -34,8 +34,8 @@ template<typename MatrixType> void submatrices(const MatrixType& m)
|
||||
*/
|
||||
|
||||
typedef typename MatrixType::Scalar Scalar;
|
||||
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
|
||||
typedef Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, 1> VectorType;
|
||||
typedef Matrix<Scalar, 1, MatrixType::Traits::ColsAtCompileTime> RowVectorType;
|
||||
int rows = m.rows();
|
||||
int cols = m.cols();
|
||||
|
||||
@ -43,9 +43,9 @@ template<typename MatrixType> void submatrices(const MatrixType& m)
|
||||
m2 = MatrixType::random(rows, cols),
|
||||
m3(rows, cols),
|
||||
mzero = MatrixType::zero(rows, cols),
|
||||
identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
identity = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::identity(rows),
|
||||
square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
|
||||
square = Matrix<Scalar, MatrixType::Traits::RowsAtCompileTime, MatrixType::Traits::RowsAtCompileTime>
|
||||
::random(rows, rows);
|
||||
VectorType v1 = VectorType::random(rows),
|
||||
v2 = VectorType::random(rows),
|
||||
|
Loading…
x
Reference in New Issue
Block a user