mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-05-22 12:37:35 +08:00
fix 2849
This commit is contained in:
parent
db85838ee2
commit
171bd08ca9
@ -603,9 +603,8 @@ class VectorwiseOp {
|
|||||||
/** Returns the expression where each subvector is the product of the vector \a other
|
/** Returns the expression where each subvector is the product of the vector \a other
|
||||||
* by the corresponding subvector of \c *this */
|
* by the corresponding subvector of \c *this */
|
||||||
template <typename OtherDerived>
|
template <typename OtherDerived>
|
||||||
EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC
|
EIGEN_DEVICE_FUNC CwiseBinaryOp<internal::scalar_product_op<Scalar, typename OtherDerived::Scalar>,
|
||||||
CwiseBinaryOp<internal::scalar_product_op<Scalar>, const ExpressionTypeNestedCleaned,
|
const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type>
|
||||||
const typename ExtendedType<OtherDerived>::Type> EIGEN_DEVICE_FUNC
|
|
||||||
operator*(const DenseBase<OtherDerived>& other) const {
|
operator*(const DenseBase<OtherDerived>& other) const {
|
||||||
EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
|
EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
|
||||||
EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType)
|
EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType)
|
||||||
@ -616,8 +615,8 @@ class VectorwiseOp {
|
|||||||
/** Returns the expression where each subvector is the quotient of the corresponding
|
/** Returns the expression where each subvector is the quotient of the corresponding
|
||||||
* subvector of \c *this by the vector \a other */
|
* subvector of \c *this by the vector \a other */
|
||||||
template <typename OtherDerived>
|
template <typename OtherDerived>
|
||||||
EIGEN_DEVICE_FUNC CwiseBinaryOp<internal::scalar_quotient_op<Scalar>, const ExpressionTypeNestedCleaned,
|
EIGEN_DEVICE_FUNC CwiseBinaryOp<internal::scalar_quotient_op<Scalar, typename OtherDerived::Scalar>,
|
||||||
const typename ExtendedType<OtherDerived>::Type>
|
const ExpressionTypeNestedCleaned, const typename ExtendedType<OtherDerived>::Type>
|
||||||
operator/(const DenseBase<OtherDerived>& other) const {
|
operator/(const DenseBase<OtherDerived>& other) const {
|
||||||
EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
|
EIGEN_STATIC_ASSERT_VECTOR_ONLY(OtherDerived)
|
||||||
EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType)
|
EIGEN_STATIC_ASSERT_ARRAYXPR(ExpressionType)
|
||||||
|
@ -214,6 +214,17 @@ void vectorwiseop_matrix(const MatrixType& m) {
|
|||||||
VERIFY_IS_EQUAL(m1.real().middleCols(0, fix<0>).colwise().maxCoeff().eval().cols(), 0);
|
VERIFY_IS_EQUAL(m1.real().middleCols(0, fix<0>).colwise().maxCoeff().eval().cols(), 0);
|
||||||
}
|
}
|
||||||
|
|
||||||
|
void vectorwiseop_mixedscalar() {
|
||||||
|
Matrix4cd a = Matrix4cd::Random();
|
||||||
|
Vector4cd b = Vector4cd::Random();
|
||||||
|
b.imag().setZero();
|
||||||
|
Vector4d b_real = b.real();
|
||||||
|
|
||||||
|
Matrix4cd c = a.array().rowwise() * b.array().transpose();
|
||||||
|
Matrix4cd d = a.array().rowwise() * b_real.array().transpose();
|
||||||
|
VERIFY_IS_CWISE_EQUAL(c, d);
|
||||||
|
}
|
||||||
|
|
||||||
EIGEN_DECLARE_TEST(vectorwiseop) {
|
EIGEN_DECLARE_TEST(vectorwiseop) {
|
||||||
CALL_SUBTEST_1(vectorwiseop_array(Array22cd()));
|
CALL_SUBTEST_1(vectorwiseop_array(Array22cd()));
|
||||||
CALL_SUBTEST_2(vectorwiseop_array(Array<double, 3, 2>()));
|
CALL_SUBTEST_2(vectorwiseop_array(Array<double, 3, 2>()));
|
||||||
@ -226,4 +237,5 @@ EIGEN_DECLARE_TEST(vectorwiseop) {
|
|||||||
MatrixXd(internal::random<int>(1, EIGEN_TEST_MAX_SIZE), internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
|
MatrixXd(internal::random<int>(1, EIGEN_TEST_MAX_SIZE), internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
|
||||||
CALL_SUBTEST_7(vectorwiseop_matrix(VectorXd(internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
|
CALL_SUBTEST_7(vectorwiseop_matrix(VectorXd(internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
|
||||||
CALL_SUBTEST_7(vectorwiseop_matrix(RowVectorXd(internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
|
CALL_SUBTEST_7(vectorwiseop_matrix(RowVectorXd(internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
|
||||||
|
CALL_SUBTEST_8(vectorwiseop_mixedscalar());
|
||||||
}
|
}
|
||||||
|
Loading…
x
Reference in New Issue
Block a user