redux.cpp 8.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008 Benoit Jacob <jacob.benoit.1@gmail.com>
  5. // Copyright (C) 2015 Gael Guennebaud <gael.guennebaud@inria.fr>
  6. //
  7. // This Source Code Form is subject to the terms of the Mozilla
  8. // Public License v. 2.0. If a copy of the MPL was not distributed
  9. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  10. #define TEST_ENABLE_TEMPORARY_TRACKING
  11. #define EIGEN_CACHEFRIENDLY_PRODUCT_THRESHOLD 8
  12. // ^^ see bug 1449
  13. #include "main.h"
  14. template<typename MatrixType> void matrixRedux(const MatrixType& m)
  15. {
  16. typedef typename MatrixType::Scalar Scalar;
  17. typedef typename MatrixType::RealScalar RealScalar;
  18. Index rows = m.rows();
  19. Index cols = m.cols();
  20. MatrixType m1 = MatrixType::Random(rows, cols);
  21. // The entries of m1 are uniformly distributed in [0,1], so m1.prod() is very small. This may lead to test
  22. // failures if we underflow into denormals. Thus, we scale so that entries are close to 1.
  23. MatrixType m1_for_prod = MatrixType::Ones(rows, cols) + RealScalar(0.2) * m1;
  24. Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> m2(rows,rows);
  25. m2.setRandom();
  26. VERIFY_IS_MUCH_SMALLER_THAN(MatrixType::Zero(rows, cols).sum(), Scalar(1));
  27. VERIFY_IS_APPROX(MatrixType::Ones(rows, cols).sum(), Scalar(float(rows*cols))); // the float() here to shut up excessive MSVC warning about int->complex conversion being lossy
  28. Scalar s(0), p(1), minc(numext::real(m1.coeff(0))), maxc(numext::real(m1.coeff(0)));
  29. for(int j = 0; j < cols; j++)
  30. for(int i = 0; i < rows; i++)
  31. {
  32. s += m1(i,j);
  33. p *= m1_for_prod(i,j);
  34. minc = (std::min)(numext::real(minc), numext::real(m1(i,j)));
  35. maxc = (std::max)(numext::real(maxc), numext::real(m1(i,j)));
  36. }
  37. const Scalar mean = s/Scalar(RealScalar(rows*cols));
  38. VERIFY_IS_APPROX(m1.sum(), s);
  39. VERIFY_IS_APPROX(m1.mean(), mean);
  40. VERIFY_IS_APPROX(m1_for_prod.prod(), p);
  41. VERIFY_IS_APPROX(m1.real().minCoeff(), numext::real(minc));
  42. VERIFY_IS_APPROX(m1.real().maxCoeff(), numext::real(maxc));
  43. // test that partial reduction works if nested expressions is forced to evaluate early
  44. VERIFY_IS_APPROX((m1.matrix() * m1.matrix().transpose()) .cwiseProduct(m2.matrix()).rowwise().sum().sum(),
  45. (m1.matrix() * m1.matrix().transpose()).eval().cwiseProduct(m2.matrix()).rowwise().sum().sum());
  46. // test slice vectorization assuming assign is ok
  47. Index r0 = internal::random<Index>(0,rows-1);
  48. Index c0 = internal::random<Index>(0,cols-1);
  49. Index r1 = internal::random<Index>(r0+1,rows)-r0;
  50. Index c1 = internal::random<Index>(c0+1,cols)-c0;
  51. VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).sum(), m1.block(r0,c0,r1,c1).eval().sum());
  52. VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).mean(), m1.block(r0,c0,r1,c1).eval().mean());
  53. VERIFY_IS_APPROX(m1_for_prod.block(r0,c0,r1,c1).prod(), m1_for_prod.block(r0,c0,r1,c1).eval().prod());
  54. VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).real().minCoeff(), m1.block(r0,c0,r1,c1).real().eval().minCoeff());
  55. VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).real().maxCoeff(), m1.block(r0,c0,r1,c1).real().eval().maxCoeff());
  56. // regression for bug 1090
  57. const int R1 = MatrixType::RowsAtCompileTime>=2 ? MatrixType::RowsAtCompileTime/2 : 6;
  58. const int C1 = MatrixType::ColsAtCompileTime>=2 ? MatrixType::ColsAtCompileTime/2 : 6;
  59. if(R1<=rows-r0 && C1<=cols-c0)
  60. {
  61. VERIFY_IS_APPROX( (m1.template block<R1,C1>(r0,c0).sum()), m1.block(r0,c0,R1,C1).sum() );
  62. }
  63. // test empty objects
  64. VERIFY_IS_APPROX(m1.block(r0,c0,0,0).sum(), Scalar(0));
  65. VERIFY_IS_APPROX(m1.block(r0,c0,0,0).prod(), Scalar(1));
  66. // test nesting complex expression
  67. VERIFY_EVALUATION_COUNT( (m1.matrix()*m1.matrix().transpose()).sum(), (MatrixType::IsVectorAtCompileTime && MatrixType::SizeAtCompileTime!=1 ? 0 : 1) );
  68. VERIFY_EVALUATION_COUNT( ((m1.matrix()*m1.matrix().transpose())+m2).sum(),(MatrixType::IsVectorAtCompileTime && MatrixType::SizeAtCompileTime!=1 ? 0 : 1));
  69. }
  70. template<typename VectorType> void vectorRedux(const VectorType& w)
  71. {
  72. using std::abs;
  73. typedef typename VectorType::Scalar Scalar;
  74. typedef typename NumTraits<Scalar>::Real RealScalar;
  75. Index size = w.size();
  76. VectorType v = VectorType::Random(size);
  77. VectorType v_for_prod = VectorType::Ones(size) + Scalar(0.2) * v; // see comment above declaration of m1_for_prod
  78. for(int i = 1; i < size; i++)
  79. {
  80. Scalar s(0), p(1);
  81. RealScalar minc(numext::real(v.coeff(0))), maxc(numext::real(v.coeff(0)));
  82. for(int j = 0; j < i; j++)
  83. {
  84. s += v[j];
  85. p *= v_for_prod[j];
  86. minc = (std::min)(minc, numext::real(v[j]));
  87. maxc = (std::max)(maxc, numext::real(v[j]));
  88. }
  89. VERIFY_IS_MUCH_SMALLER_THAN(abs(s - v.head(i).sum()), Scalar(1));
  90. VERIFY_IS_APPROX(p, v_for_prod.head(i).prod());
  91. VERIFY_IS_APPROX(minc, v.real().head(i).minCoeff());
  92. VERIFY_IS_APPROX(maxc, v.real().head(i).maxCoeff());
  93. }
  94. for(int i = 0; i < size-1; i++)
  95. {
  96. Scalar s(0), p(1);
  97. RealScalar minc(numext::real(v.coeff(i))), maxc(numext::real(v.coeff(i)));
  98. for(int j = i; j < size; j++)
  99. {
  100. s += v[j];
  101. p *= v_for_prod[j];
  102. minc = (std::min)(minc, numext::real(v[j]));
  103. maxc = (std::max)(maxc, numext::real(v[j]));
  104. }
  105. VERIFY_IS_MUCH_SMALLER_THAN(abs(s - v.tail(size-i).sum()), Scalar(1));
  106. VERIFY_IS_APPROX(p, v_for_prod.tail(size-i).prod());
  107. VERIFY_IS_APPROX(minc, v.real().tail(size-i).minCoeff());
  108. VERIFY_IS_APPROX(maxc, v.real().tail(size-i).maxCoeff());
  109. }
  110. for(int i = 0; i < size/2; i++)
  111. {
  112. Scalar s(0), p(1);
  113. RealScalar minc(numext::real(v.coeff(i))), maxc(numext::real(v.coeff(i)));
  114. for(int j = i; j < size-i; j++)
  115. {
  116. s += v[j];
  117. p *= v_for_prod[j];
  118. minc = (std::min)(minc, numext::real(v[j]));
  119. maxc = (std::max)(maxc, numext::real(v[j]));
  120. }
  121. VERIFY_IS_MUCH_SMALLER_THAN(abs(s - v.segment(i, size-2*i).sum()), Scalar(1));
  122. VERIFY_IS_APPROX(p, v_for_prod.segment(i, size-2*i).prod());
  123. VERIFY_IS_APPROX(minc, v.real().segment(i, size-2*i).minCoeff());
  124. VERIFY_IS_APPROX(maxc, v.real().segment(i, size-2*i).maxCoeff());
  125. }
  126. // test empty objects
  127. VERIFY_IS_APPROX(v.head(0).sum(), Scalar(0));
  128. VERIFY_IS_APPROX(v.tail(0).prod(), Scalar(1));
  129. VERIFY_RAISES_ASSERT(v.head(0).mean());
  130. VERIFY_RAISES_ASSERT(v.head(0).minCoeff());
  131. VERIFY_RAISES_ASSERT(v.head(0).maxCoeff());
  132. }
  133. EIGEN_DECLARE_TEST(redux)
  134. {
  135. // the max size cannot be too large, otherwise reduxion operations obviously generate large errors.
  136. int maxsize = (std::min)(100,EIGEN_TEST_MAX_SIZE);
  137. TEST_SET_BUT_UNUSED_VARIABLE(maxsize);
  138. for(int i = 0; i < g_repeat; i++) {
  139. CALL_SUBTEST_1( matrixRedux(Matrix<float, 1, 1>()) );
  140. CALL_SUBTEST_1( matrixRedux(Array<float, 1, 1>()) );
  141. CALL_SUBTEST_2( matrixRedux(Matrix2f()) );
  142. CALL_SUBTEST_2( matrixRedux(Array2f()) );
  143. CALL_SUBTEST_2( matrixRedux(Array22f()) );
  144. CALL_SUBTEST_3( matrixRedux(Matrix4d()) );
  145. CALL_SUBTEST_3( matrixRedux(Array4d()) );
  146. CALL_SUBTEST_3( matrixRedux(Array44d()) );
  147. CALL_SUBTEST_4( matrixRedux(MatrixXcf(internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  148. CALL_SUBTEST_4( matrixRedux(ArrayXXcf(internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  149. CALL_SUBTEST_5( matrixRedux(MatrixXd (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  150. CALL_SUBTEST_5( matrixRedux(ArrayXXd (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  151. CALL_SUBTEST_6( matrixRedux(MatrixXi (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  152. CALL_SUBTEST_6( matrixRedux(ArrayXXi (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  153. }
  154. for(int i = 0; i < g_repeat; i++) {
  155. CALL_SUBTEST_7( vectorRedux(Vector4f()) );
  156. CALL_SUBTEST_7( vectorRedux(Array4f()) );
  157. CALL_SUBTEST_5( vectorRedux(VectorXd(internal::random<int>(1,maxsize))) );
  158. CALL_SUBTEST_5( vectorRedux(ArrayXd(internal::random<int>(1,maxsize))) );
  159. CALL_SUBTEST_8( vectorRedux(VectorXf(internal::random<int>(1,maxsize))) );
  160. CALL_SUBTEST_8( vectorRedux(ArrayXf(internal::random<int>(1,maxsize))) );
  161. }
  162. }