- make RowsAtCompileTime and ColsAtCompileTime public in

MatrixBase and private in derived types
- initial documentation in MatrixBase
This commit is contained in:
Benoit Jacob 2007-12-19 09:30:53 +00:00
parent 59be5c3124
commit cddeeee17d
23 changed files with 157 additions and 96 deletions

View File

@ -35,9 +35,6 @@ template<typename MatrixType, int BlockRows, int BlockCols> class Block
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, Block<MatrixType, BlockRows, BlockCols> >;
static const int RowsAtCompileTime = BlockRows,
ColsAtCompileTime = BlockCols;
Block(const MatRef& matrix, int startRow, int startCol)
: m_matrix(matrix), m_startRow(startRow), m_startCol(startCol)
{
@ -52,6 +49,9 @@ template<typename MatrixType, int BlockRows, int BlockCols> class Block
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Block)
private:
static const int _RowsAtCompileTime = BlockRows,
_ColsAtCompileTime = BlockCols;
const Block& _ref() const { return *this; }
int _rows() const { return BlockRows; }
int _cols() const { return BlockCols; }

View File

@ -34,15 +34,14 @@ template<typename NewScalar, typename MatrixType> class Cast : NoOperatorEquals,
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, Cast<Scalar, MatrixType> >;
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
Cast(const MatRef& matrix) : m_matrix(matrix) {}
Cast(const Cast& other)
: m_matrix(other.m_matrix) {}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Cast& _ref() const { return *this; }
int _rows() const { return m_matrix.rows(); }
int _cols() const { return m_matrix.cols(); }

View File

@ -34,9 +34,6 @@ template<typename MatrixType> class Column
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, Column<MatrixType> >;
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = 1;
Column(const MatRef& matrix, int col)
: m_matrix(matrix), m_col(col)
{
@ -49,6 +46,8 @@ template<typename MatrixType> class Column
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Column)
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = 1;
const Column& _ref() const { return *this; }
int _rows() const { return m_matrix.rows(); }
int _cols() const { return 1; }

View File

@ -34,15 +34,15 @@ template<typename MatrixType> class Conjugate : NoOperatorEquals,
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, Conjugate<MatrixType> >;
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
Conjugate(const MatRef& matrix) : m_matrix(matrix) {}
Conjugate(const Conjugate& other)
: m_matrix(other.m_matrix) {}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Conjugate& _ref() const { return *this; }
int _rows() const { return m_matrix.rows(); }
int _cols() const { return m_matrix.cols(); }

View File

@ -34,9 +34,6 @@ template<typename MatrixType> class DiagonalCoeffs
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, DiagonalCoeffs<MatrixType> >;
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = 1;
DiagonalCoeffs(const MatRef& matrix) : m_matrix(matrix) {}
DiagonalCoeffs(const DiagonalCoeffs& other) : m_matrix(other.m_matrix) {}
@ -44,6 +41,9 @@ template<typename MatrixType> class DiagonalCoeffs
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(DiagonalCoeffs)
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = 1;
const DiagonalCoeffs& _ref() const { return *this; }
int _rows() const { return std::min(m_matrix.rows(), m_matrix.cols()); }
int _cols() const { return 1; }

View File

@ -36,18 +36,18 @@ class DiagonalMatrix : NoOperatorEquals,
typedef typename CoeffsVectorType::Ref CoeffsVecRef;
friend class MatrixBase<Scalar, DiagonalMatrix<MatrixType, CoeffsVectorType> >;
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
DiagonalMatrix(const CoeffsVecRef& coeffs) : m_coeffs(coeffs)
{
assert(CoeffsVectorType::IsVector
&& RowsAtCompileTime == ColsAtCompileTime
&& RowsAtCompileTime == CoeffsVectorType::SizeAtCompileTime
&& _RowsAtCompileTime == _ColsAtCompileTime
&& _RowsAtCompileTime == CoeffsVectorType::SizeAtCompileTime
&& coeffs.size() > 0);
}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const DiagonalMatrix& _ref() const { return *this; }
int _rows() const { return m_coeffs.size(); }
int _cols() const { return m_coeffs.size(); }

View File

@ -35,9 +35,6 @@ template<typename Lhs, typename Rhs> class Difference : NoOperatorEquals,
typedef typename Rhs::Ref RhsRef;
friend class MatrixBase<Scalar, Difference>;
static const int RowsAtCompileTime = Lhs::RowsAtCompileTime,
ColsAtCompileTime = Rhs::ColsAtCompileTime;
Difference(const LhsRef& lhs, const RhsRef& rhs)
: m_lhs(lhs), m_rhs(rhs)
{
@ -48,6 +45,9 @@ template<typename Lhs, typename Rhs> class Difference : NoOperatorEquals,
: m_lhs(other.m_lhs), m_rhs(other.m_rhs) {}
private:
static const int _RowsAtCompileTime = Lhs::RowsAtCompileTime,
_ColsAtCompileTime = Rhs::ColsAtCompileTime;
const Difference& _ref() const { return *this; }
int _rows() const { return m_lhs.rows(); }
int _cols() const { return m_lhs.cols(); }

View File

@ -34,10 +34,6 @@ template<typename MatrixType> class DynBlock
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, DynBlock<MatrixType> >;
static const int
RowsAtCompileTime = MatrixType::RowsAtCompileTime == 1 ? 1 : Dynamic,
ColsAtCompileTime = MatrixType::ColsAtCompileTime == 1 ? 1 : Dynamic;
DynBlock(const MatRef& matrix,
int startRow, int startCol,
int blockRows, int blockCols)
@ -56,6 +52,10 @@ template<typename MatrixType> class DynBlock
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(DynBlock)
private:
static const int
_RowsAtCompileTime = MatrixType::RowsAtCompileTime == 1 ? 1 : Dynamic,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime == 1 ? 1 : Dynamic;
const DynBlock& _ref() const { return *this; }
int _rows() const { return m_blockRows; }
int _cols() const { return m_blockCols; }

View File

@ -33,15 +33,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);
assert(rows > 0 && _RowsAtCompileTime == _ColsAtCompileTime);
}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Identity& _ref() const { return *this; }
int _rows() const { return m_rows; }
int _cols() const { return m_rows; }

View File

@ -23,8 +23,8 @@
// License. This exception does not invalidate any other reasons why a work
// based on this file might be covered by the GNU General Public License.
#ifndef EIGEN_FROMARRAY_H
#define EIGEN_FROMARRAY_H
#ifndef EIGEN_MAP_H
#define EIGEN_MAP_H
template<typename MatrixType> class Map
: public MatrixBase<typename MatrixType::Scalar, Map<MatrixType> >
@ -33,9 +33,6 @@ 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;
Map(const Scalar* data, int rows, int cols) : m_data(data), m_rows(rows), m_cols(cols)
{
assert(rows > 0 && cols > 0);
@ -44,6 +41,9 @@ template<typename MatrixType> class Map
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Map)
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Map& _ref() const { return *this; }
int _rows() const { return m_rows; }
int _cols() const { return m_cols; }
@ -131,4 +131,4 @@ Matrix<_Scalar, _Rows, _Cols>
*this = map(data);
}
#endif // EIGEN_FROMARRAY_H
#endif // EIGEN_MAP_H

View File

@ -38,8 +38,6 @@ class Matrix : public MatrixBase<_Scalar, Matrix<_Scalar, _Rows, _Cols> >,
typedef MatrixRef<Matrix> Ref;
friend class MatrixRef<Matrix>;
static const int RowsAtCompileTime = _Rows, ColsAtCompileTime = _Cols;
const Scalar* data() const
{ return Storage::m_data; }
@ -47,6 +45,8 @@ class Matrix : public MatrixBase<_Scalar, Matrix<_Scalar, _Rows, _Cols> >,
{ return Storage::m_data; }
private:
static const int _RowsAtCompileTime = _Rows, _ColsAtCompileTime = _Cols;
Ref _ref() const { return Ref(*this); }
const Scalar& _coeff(int row, int col) const
@ -80,15 +80,15 @@ class Matrix : public MatrixBase<_Scalar, Matrix<_Scalar, _Rows, _Cols> >,
explicit Matrix() : Storage()
{
assert(RowsAtCompileTime > 0 && ColsAtCompileTime > 0);
assert(_RowsAtCompileTime > 0 && _ColsAtCompileTime > 0);
}
explicit Matrix(int dim) : Storage(dim)
{
assert(dim > 0);
assert((RowsAtCompileTime == 1
&& (ColsAtCompileTime == Dynamic || ColsAtCompileTime == dim))
|| (ColsAtCompileTime == 1
&& (RowsAtCompileTime == Dynamic || RowsAtCompileTime == dim)));
assert((_RowsAtCompileTime == 1
&& (_ColsAtCompileTime == Dynamic || _ColsAtCompileTime == dim))
|| (_ColsAtCompileTime == 1
&& (_RowsAtCompileTime == Dynamic || _RowsAtCompileTime == dim)));
}
// this constructor is very tricky.
@ -101,44 +101,44 @@ class Matrix : public MatrixBase<_Scalar, Matrix<_Scalar, _Rows, _Cols> >,
// does what we want to, so it only remains to add some asserts.
Matrix(int x, int y) : Storage(x, y)
{
if((RowsAtCompileTime == 1 && ColsAtCompileTime == 2)
|| (RowsAtCompileTime == 2 && ColsAtCompileTime == 1))
if((_RowsAtCompileTime == 1 && _ColsAtCompileTime == 2)
|| (_RowsAtCompileTime == 2 && _ColsAtCompileTime == 1))
{
(Storage::m_data)[0] = x;
(Storage::m_data)[1] = y;
}
else
{
assert(x > 0 && (RowsAtCompileTime == Dynamic || RowsAtCompileTime == x)
&& y > 0 && (ColsAtCompileTime == Dynamic || ColsAtCompileTime == y));
assert(x > 0 && (_RowsAtCompileTime == Dynamic || _RowsAtCompileTime == x)
&& y > 0 && (_ColsAtCompileTime == Dynamic || _ColsAtCompileTime == y));
}
}
Matrix(const float& x, const float& y)
{
assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 2)
|| (RowsAtCompileTime == 2 && ColsAtCompileTime == 1));
assert((_RowsAtCompileTime == 1 && _ColsAtCompileTime == 2)
|| (_RowsAtCompileTime == 2 && _ColsAtCompileTime == 1));
(Storage::m_data)[0] = x;
(Storage::m_data)[1] = y;
}
Matrix(const double& x, const double& y)
{
assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 2)
|| (RowsAtCompileTime == 2 && ColsAtCompileTime == 1));
assert((_RowsAtCompileTime == 1 && _ColsAtCompileTime == 2)
|| (_RowsAtCompileTime == 2 && _ColsAtCompileTime == 1));
(Storage::m_data)[0] = x;
(Storage::m_data)[1] = y;
}
Matrix(const Scalar& x, const Scalar& y, const Scalar& z)
{
assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 3)
|| (RowsAtCompileTime == 3 && ColsAtCompileTime == 1));
assert((_RowsAtCompileTime == 1 && _ColsAtCompileTime == 3)
|| (_RowsAtCompileTime == 3 && _ColsAtCompileTime == 1));
(Storage::m_data)[0] = x;
(Storage::m_data)[1] = y;
(Storage::m_data)[2] = z;
}
Matrix(const Scalar& x, const Scalar& y, const Scalar& z, const Scalar& w)
{
assert((RowsAtCompileTime == 1 && ColsAtCompileTime == 4)
|| (RowsAtCompileTime == 4 && ColsAtCompileTime == 1));
assert((_RowsAtCompileTime == 1 && _ColsAtCompileTime == 4)
|| (_RowsAtCompileTime == 4 && _ColsAtCompileTime == 1));
(Storage::m_data)[0] = x;
(Storage::m_data)[1] = y;
(Storage::m_data)[2] = z;

View File

@ -26,41 +26,104 @@
#ifndef EIGEN_MATRIXBASE_H
#define EIGEN_MATRIXBASE_H
/** \class MatrixBase
*
* \brief Base class for all matrices, vectors, and expressions
*
* This class is the base that is inherited by all matrix, vector, and expression
* types. Most of the Eigen API is contained in this class.
*
* This class takes two template parameters:
* \param Scalar the type of the coefficients, e.g. float, double, etc.
* \param Derived the derived type, e.g. a matrix type, or an expression, etc.
* Indeed, a separate MatrixBase type is generated for each derived type
* so one knows from inside MatrixBase, at compile-time, what the derived type is.
*
* When writing a function taking Eigen objects as argument, if you want your function
* to take as argument any matrix, vector, or expression, just let it take a
* MatrixBase argument. As an example, here is a function printFirstRow which, given
* a matrix, vector, or expression \a x, prints the first row of \a x.
*
* \code
template<typename Scalar, typename Derived>
void printFirstRow(const Eigen::MatrixBase<Scalar, Derived>& m)
{
cout << x.row(0) << endl;
}
* \endcode
*/
template<typename Scalar, typename Derived> class MatrixBase
{
static const int RowsAtCompileTime = Derived::RowsAtCompileTime,
ColsAtCompileTime = Derived::ColsAtCompileTime;
public:
/** The number of rows and of columns at compile-time. These are just
* copies of the values 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(), SizeAtCompileTime */
static const int RowsAtCompileTime = Derived::_RowsAtCompileTime,
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
= RowsAtCompileTime == Dynamic || ColsAtCompileTime == Dynamic
? Dynamic : RowsAtCompileTime * 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 IsVector = RowsAtCompileTime == 1 || 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
* is just this MatrixBase type itself, i.e. expressions are just passed by value
* and the compiler is supposed to be clever enough to optimize that. If, on the
* other hand, this MatrixBase type is an actual matrix or vector, then \a Ref is
* a typedef MatrixRef, which is like a reference, so that matrices and vectors
* are passed by reference, not by value. \sa ref()*/
typedef typename ForwardDecl<Derived>::Ref Ref;
/** This is the "real scalar" type; if the \a Scalar type is already real numbers
* (e.g. int, float or double) then RealScalar is just the same as \a Scalar. If
* \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 */
int rows() const { return static_cast<const Derived *>(this)->_rows(); }
/** \returns the number of columns. \sa row(), ColsAtCompileTime*/
int cols() const { return static_cast<const Derived *>(this)->_cols(); }
/** \returns the number of coefficients, which is \a rows()*cols().
* \sa rows(), cols(). */
int size() const { return rows() * cols(); }
/** \returns a Ref to *this. \sa Ref */
Ref ref() const
{ return static_cast<const Derived *>(this)->_ref(); }
/** Copies \a other into *this. \returns a reference to *this. */
template<typename OtherDerived>
Derived& operator=(const MatrixBase<Scalar, OtherDerived>& other);
//special case of the above template operator=, in order to prevent the compiler
// Special case of the above template operator=, in order to prevent the compiler
//from generating a default operator= (issue hit with g++ 4.1)
Derived& operator=(const MatrixBase& other)
{
return this->operator=<Derived>(other);
}
/** \returns an expression of *this with the \a Scalar type casted to
* \a NewScalar. */
template<typename NewScalar> const Cast<NewScalar, Derived> cast() const;
/** \returns an expression of the \a i-th row of *this.
* \sa col(int)*/
Row<Derived> row(int i) const;
/** \returns an expression of the \a i-th column of *this.
* \sa row(int) */
Column<Derived> col(int i) const;
/** \return an expression of the (\a row, \a col)-minor of *this,
* i.e. an expression constructed from *this by removing the specified
* row and column. */
Minor<Derived> minor(int row, int col) const;
DynBlock<Derived> dynBlock(int startRow, int startCol,

View File

@ -33,9 +33,6 @@ template<typename MatrixType> class MatrixRef
typedef typename MatrixType::Scalar Scalar;
friend class MatrixBase<Scalar, MatrixRef>;
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
MatrixRef(const MatrixType& matrix) : m_matrix(*const_cast<MatrixType*>(&matrix)) {}
MatrixRef(const MatrixRef& other) : m_matrix(other.m_matrix) {}
~MatrixRef() {}
@ -43,6 +40,9 @@ template<typename MatrixType> class MatrixRef
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(MatrixRef)
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
int _rows() const { return m_matrix.rows(); }
int _cols() const { return m_matrix.cols(); }

View File

@ -34,12 +34,6 @@ template<typename MatrixType> class Minor
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, Minor<MatrixType> >;
static const int
RowsAtCompileTime = (MatrixType::RowsAtCompileTime != Dynamic) ?
MatrixType::RowsAtCompileTime - 1 : Dynamic,
ColsAtCompileTime = (MatrixType::ColsAtCompileTime != Dynamic) ?
MatrixType::ColsAtCompileTime - 1 : Dynamic;
Minor(const MatRef& matrix,
int row, int col)
: m_matrix(matrix), m_row(row), m_col(col)
@ -54,6 +48,12 @@ template<typename MatrixType> class Minor
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Minor)
private:
static const int
_RowsAtCompileTime = (MatrixType::RowsAtCompileTime != Dynamic) ?
MatrixType::RowsAtCompileTime - 1 : Dynamic,
_ColsAtCompileTime = (MatrixType::ColsAtCompileTime != Dynamic) ?
MatrixType::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; }

View File

@ -33,15 +33,15 @@ 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;
Ones(int rows, int cols) : m_rows(rows), m_cols(cols)
{
assert(rows > 0 && cols > 0);
}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Ones& _ref() const { return *this; }
int _rows() const { return m_rows; }
int _cols() const { return m_cols; }

View File

@ -34,15 +34,15 @@ template<typename MatrixType> class Opposite : NoOperatorEquals,
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, Opposite<MatrixType> >;
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
Opposite(const MatRef& matrix) : m_matrix(matrix) {}
Opposite(const Opposite& other)
: m_matrix(other.m_matrix) {}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Opposite& _ref() const { return *this; }
int _rows() const { return m_matrix.rows(); }
int _cols() const { return m_matrix.cols(); }

View File

@ -69,9 +69,6 @@ template<typename Lhs, typename Rhs> class Product : NoOperatorEquals,
typedef typename Rhs::Ref RhsRef;
friend class MatrixBase<Scalar, Product>;
static const int RowsAtCompileTime = Lhs::RowsAtCompileTime,
ColsAtCompileTime = Rhs::ColsAtCompileTime;
Product(const LhsRef& lhs, const RhsRef& rhs)
: m_lhs(lhs), m_rhs(rhs)
{
@ -82,6 +79,9 @@ template<typename Lhs, typename Rhs> class Product : NoOperatorEquals,
: m_lhs(other.m_lhs), m_rhs(other.m_rhs) {}
private:
static const int _RowsAtCompileTime = Lhs::RowsAtCompileTime,
_ColsAtCompileTime = Rhs::ColsAtCompileTime;
const Product& _ref() const { return *this; }
int _rows() const { return m_lhs.rows(); }
int _cols() const { return m_rhs.cols(); }

View File

@ -33,15 +33,15 @@ 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;
Random(int rows, int cols) : m_rows(rows), m_cols(cols)
{
assert(rows > 0 && cols > 0);
}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Random& _ref() const { return *this; }
int _rows() const { return m_rows; }
int _cols() const { return m_cols; }

View File

@ -34,9 +34,6 @@ template<typename MatrixType> class Row
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, Row<MatrixType> >;
static const int RowsAtCompileTime = 1,
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
Row(const MatRef& matrix, int row)
: m_matrix(matrix), m_row(row)
{
@ -55,6 +52,9 @@ template<typename MatrixType> class Row
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Row)
private:
static const int _RowsAtCompileTime = 1,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Row& _ref() const { return *this; }
int _rows() const { return 1; }

View File

@ -34,9 +34,6 @@ template<typename MatrixType> class ScalarMultiple : NoOperatorEquals,
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<typename MatrixType::Scalar, ScalarMultiple<MatrixType> >;
static const int RowsAtCompileTime = MatrixType::RowsAtCompileTime,
ColsAtCompileTime = MatrixType::ColsAtCompileTime;
ScalarMultiple(const MatRef& matrix, Scalar scalar)
: m_matrix(matrix), m_scalar(scalar) {}
@ -44,6 +41,9 @@ template<typename MatrixType> class ScalarMultiple : NoOperatorEquals,
: m_matrix(other.m_matrix), m_scalar(other.m_scalar) {}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const ScalarMultiple& _ref() const { return *this; }
int _rows() const { return m_matrix.rows(); }
int _cols() const { return m_matrix.cols(); }

View File

@ -35,9 +35,6 @@ template<typename Lhs, typename Rhs> class Sum : NoOperatorEquals,
typedef typename Rhs::Ref RhsRef;
friend class MatrixBase<Scalar, Sum>;
static const int RowsAtCompileTime = Lhs::RowsAtCompileTime,
ColsAtCompileTime = Rhs::ColsAtCompileTime;
Sum(const LhsRef& lhs, const RhsRef& rhs)
: m_lhs(lhs), m_rhs(rhs)
{
@ -47,6 +44,9 @@ template<typename Lhs, typename Rhs> class Sum : NoOperatorEquals,
Sum(const Sum& other) : m_lhs(other.m_lhs), m_rhs(other.m_rhs) {}
private:
static const int _RowsAtCompileTime = Lhs::RowsAtCompileTime,
_ColsAtCompileTime = Rhs::ColsAtCompileTime;
const Sum& _ref() const { return *this; }
int _rows() const { return m_lhs.rows(); }
int _cols() const { return m_lhs.cols(); }

View File

@ -34,9 +34,6 @@ template<typename MatrixType> class Transpose
typedef typename MatrixType::Ref MatRef;
friend class MatrixBase<Scalar, Transpose<MatrixType> >;
static const int RowsAtCompileTime = MatrixType::ColsAtCompileTime,
ColsAtCompileTime = MatrixType::RowsAtCompileTime;
Transpose(const MatRef& matrix) : m_matrix(matrix) {}
Transpose(const Transpose& other)
@ -45,6 +42,9 @@ template<typename MatrixType> class Transpose
EIGEN_INHERIT_ASSIGNMENT_OPERATORS(Transpose)
private:
static const int _RowsAtCompileTime = MatrixType::ColsAtCompileTime,
_ColsAtCompileTime = MatrixType::RowsAtCompileTime;
const Transpose& _ref() const { return *this; }
int _rows() const { return m_matrix.cols(); }
int _cols() const { return m_matrix.rows(); }

View File

@ -33,15 +33,15 @@ 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;
Zero(int rows, int cols) : m_rows(rows), m_cols(cols)
{
assert(rows > 0 && cols > 0);
}
private:
static const int _RowsAtCompileTime = MatrixType::RowsAtCompileTime,
_ColsAtCompileTime = MatrixType::ColsAtCompileTime;
const Zero& _ref() const { return *this; }
int _rows() const { return m_rows; }
int _cols() const { return m_cols; }