mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-21 09:09:36 +08:00
bug #1376: add missing assertion on size mismatch with compound assignment operators (e.g., mat += mat.col(j))
This commit is contained in:
parent
b0db4eff36
commit
ba3f977946
@ -701,6 +701,26 @@ protected:
|
||||
* Part 5 : Entry point for dense rectangular assignment
|
||||
***************************************************************************/
|
||||
|
||||
template<typename DstXprType,typename SrcXprType, typename Functor>
|
||||
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||
void resize_if_allowed(DstXprType &dst, const SrcXprType& src, const Functor &/*func*/)
|
||||
{
|
||||
EIGEN_ONLY_USED_FOR_DEBUG(dst);
|
||||
EIGEN_ONLY_USED_FOR_DEBUG(src);
|
||||
eigen_assert(dst.rows() == src.rows() && dst.cols() == src.cols());
|
||||
};
|
||||
|
||||
template<typename DstXprType,typename SrcXprType, typename T1, typename T2>
|
||||
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE
|
||||
void resize_if_allowed(DstXprType &dst, const SrcXprType& src, const internal::assign_op<T1,T2> &/*func*/)
|
||||
{
|
||||
Index dstRows = src.rows();
|
||||
Index dstCols = src.cols();
|
||||
if(((dst.rows()!=dstRows) || (dst.cols()!=dstCols)))
|
||||
dst.resize(dstRows, dstCols);
|
||||
eigen_assert(dst.rows() == dstRows && dst.cols() == dstCols);
|
||||
};
|
||||
|
||||
template<typename DstXprType, typename SrcXprType, typename Functor>
|
||||
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void call_dense_assignment_loop(DstXprType& dst, const SrcXprType& src, const Functor &func)
|
||||
{
|
||||
@ -711,10 +731,7 @@ EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void call_dense_assignment_loop(DstXprType
|
||||
|
||||
// NOTE To properly handle A = (A*A.transpose())/s with A rectangular,
|
||||
// we need to resize the destination after the source evaluator has been created.
|
||||
Index dstRows = src.rows();
|
||||
Index dstCols = src.cols();
|
||||
if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
|
||||
dst.resize(dstRows, dstCols);
|
||||
resize_if_allowed(dst, src, func);
|
||||
|
||||
DstEvaluatorType dstEvaluator(dst);
|
||||
|
||||
|
@ -28,7 +28,7 @@ template<typename DstScalar,typename SrcScalar> struct assign_op {
|
||||
{ internal::pstoret<DstScalar,Packet,Alignment>(a,b); }
|
||||
};
|
||||
|
||||
// Empty overload for void type (used by PermutationMatrix
|
||||
// Empty overload for void type (used by PermutationMatrix)
|
||||
template<typename DstScalar> struct assign_op<DstScalar,void> {};
|
||||
|
||||
template<typename DstScalar,typename SrcScalar>
|
||||
|
@ -143,10 +143,7 @@ struct Assignment<DstXprType, SrcXprType, Functor, Sparse2Dense>
|
||||
dst.setZero();
|
||||
|
||||
internal::evaluator<SrcXprType> srcEval(src);
|
||||
Index dstRows = src.rows();
|
||||
Index dstCols = src.cols();
|
||||
if((dst.rows()!=dstRows) || (dst.cols()!=dstCols))
|
||||
dst.resize(dstRows, dstCols);
|
||||
resize_if_allowed(dst, src, func);
|
||||
internal::evaluator<DstXprType> dstEval(dst);
|
||||
|
||||
const Index outerEvaluationSize = (internal::evaluator<SrcXprType>::Flags&RowMajorBit) ? src.rows() : src.cols();
|
||||
|
@ -186,6 +186,14 @@ template<typename MatrixType> void block(const MatrixType& m)
|
||||
VERIFY_IS_EQUAL( (m1.template block<1,Dynamic>(0,1,1,0)), m1.block(0,1,1,0));
|
||||
VERIFY_IS_EQUAL( ((m1*1).template block<Dynamic,1>(1,0,0,1)), m1.block(1,0,0,1));
|
||||
VERIFY_IS_EQUAL( ((m1*1).template block<1,Dynamic>(0,1,1,0)), m1.block(0,1,1,0));
|
||||
|
||||
if (rows>=2 && cols>=2)
|
||||
{
|
||||
VERIFY_RAISES_ASSERT( m1 += m1.col(0) );
|
||||
VERIFY_RAISES_ASSERT( m1 -= m1.col(0) );
|
||||
VERIFY_RAISES_ASSERT( m1.array() *= m1.col(0).array() );
|
||||
VERIFY_RAISES_ASSERT( m1.array() /= m1.col(0).array() );
|
||||
}
|
||||
}
|
||||
|
||||
|
||||
|
@ -68,6 +68,21 @@ template<typename MatrixType> void diagonal(const MatrixType& m)
|
||||
}
|
||||
}
|
||||
|
||||
template<typename MatrixType> void diagonal_assert(const MatrixType& m) {
|
||||
Index rows = m.rows();
|
||||
Index cols = m.cols();
|
||||
|
||||
MatrixType m1 = MatrixType::Random(rows, cols);
|
||||
|
||||
if (rows>=2 && cols>=2)
|
||||
{
|
||||
VERIFY_RAISES_ASSERT( m1 += m1.diagonal() );
|
||||
VERIFY_RAISES_ASSERT( m1 -= m1.diagonal() );
|
||||
VERIFY_RAISES_ASSERT( m1.array() *= m1.diagonal().array() );
|
||||
VERIFY_RAISES_ASSERT( m1.array() /= m1.diagonal().array() );
|
||||
}
|
||||
}
|
||||
|
||||
void test_diagonal()
|
||||
{
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
@ -81,4 +96,6 @@ void test_diagonal()
|
||||
CALL_SUBTEST_1( diagonal(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
|
||||
CALL_SUBTEST_1( diagonal(Matrix<float,Dynamic,4>(3, 4)) );
|
||||
}
|
||||
|
||||
CALL_SUBTEST_1( diagonal_assert(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
|
||||
}
|
||||
|
@ -214,6 +214,14 @@ template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& re
|
||||
VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
|
||||
VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
|
||||
|
||||
if (rows>=2 && cols>=2)
|
||||
{
|
||||
VERIFY_RAISES_ASSERT( m1 += m1.innerVector(0) );
|
||||
VERIFY_RAISES_ASSERT( m1 -= m1.innerVector(0) );
|
||||
VERIFY_RAISES_ASSERT( refM1 -= m1.innerVector(0) );
|
||||
VERIFY_RAISES_ASSERT( refM1 += m1.innerVector(0) );
|
||||
}
|
||||
|
||||
// test aliasing
|
||||
VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
|
||||
VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
|
||||
@ -428,7 +436,7 @@ template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& re
|
||||
m3 = m2.template triangularView<StrictlyLower>();
|
||||
VERIFY_IS_APPROX(m3, refMat3);
|
||||
|
||||
// check sparse-traingular to dense
|
||||
// check sparse-triangular to dense
|
||||
refMat3 = m2.template triangularView<StrictlyUpper>();
|
||||
VERIFY_IS_APPROX(refMat3, DenseMatrix(refMat2.template triangularView<StrictlyUpper>()));
|
||||
}
|
||||
|
Loading…
x
Reference in New Issue
Block a user