mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-08-11 11:19:02 +08:00
make eigen2 submatrices test pass
This commit is contained in:
parent
cc2b7a5397
commit
5be269db88
@ -56,7 +56,7 @@ struct traits<Minor<MatrixType> >
|
|||||||
int(MatrixType::MaxRowsAtCompileTime) - 1 : Dynamic,
|
int(MatrixType::MaxRowsAtCompileTime) - 1 : Dynamic,
|
||||||
MaxColsAtCompileTime = (MatrixType::MaxColsAtCompileTime != Dynamic) ?
|
MaxColsAtCompileTime = (MatrixType::MaxColsAtCompileTime != Dynamic) ?
|
||||||
int(MatrixType::MaxColsAtCompileTime) - 1 : Dynamic,
|
int(MatrixType::MaxColsAtCompileTime) - 1 : Dynamic,
|
||||||
Flags = _MatrixTypeNested::Flags & HereditaryBits,
|
Flags = _MatrixTypeNested::Flags & (HereditaryBits | LvalueBit),
|
||||||
CoeffReadCost = _MatrixTypeNested::CoeffReadCost // minor is used typically on tiny matrices,
|
CoeffReadCost = _MatrixTypeNested::CoeffReadCost // minor is used typically on tiny matrices,
|
||||||
// where loops are unrolled and the 'if' evaluates at compile time
|
// where loops are unrolled and the 'if' evaluates at compile time
|
||||||
};
|
};
|
||||||
|
@ -111,8 +111,10 @@ template<typename MatrixType> void submatrices(const MatrixType& m)
|
|||||||
m2.diagonal()[0] *= 3;
|
m2.diagonal()[0] *= 3;
|
||||||
VERIFY_IS_APPROX(m2.diagonal()[0], static_cast<Scalar>(6) * m1.diagonal()[0]);
|
VERIFY_IS_APPROX(m2.diagonal()[0], static_cast<Scalar>(6) * m1.diagonal()[0]);
|
||||||
|
|
||||||
const int BlockRows = EIGEN_ENUM_MIN(MatrixType::RowsAtCompileTime,2);
|
enum {
|
||||||
const int BlockCols = EIGEN_ENUM_MIN(MatrixType::ColsAtCompileTime,5);
|
BlockRows = EIGEN_SIZE_MIN_PREFER_FIXED(MatrixType::RowsAtCompileTime,2),
|
||||||
|
BlockCols = EIGEN_SIZE_MIN_PREFER_FIXED(MatrixType::ColsAtCompileTime,5)
|
||||||
|
};
|
||||||
if (rows>=5 && cols>=8)
|
if (rows>=5 && cols>=8)
|
||||||
{
|
{
|
||||||
// test fixed block() as lvalue
|
// test fixed block() as lvalue
|
||||||
|
Loading…
x
Reference in New Issue
Block a user