You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

159 lines
6.8 KiB

  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. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "main.h"
  10. template<typename MatrixType> void matrixRedux(const MatrixType& m)
  11. {
  12. typedef typename MatrixType::Index Index;
  13. typedef typename MatrixType::Scalar Scalar;
  14. typedef typename MatrixType::RealScalar RealScalar;
  15. Index rows = m.rows();
  16. Index cols = m.cols();
  17. MatrixType m1 = MatrixType::Random(rows, cols);
  18. // The entries of m1 are uniformly distributed in [0,1], so m1.prod() is very small. This may lead to test
  19. // failures if we underflow into denormals. Thus, we scale so that entires are close to 1.
  20. MatrixType m1_for_prod = MatrixType::Ones(rows, cols) + RealScalar(0.2) * m1;
  21. VERIFY_IS_MUCH_SMALLER_THAN(MatrixType::Zero(rows, cols).sum(), Scalar(1));
  22. 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
  23. Scalar s(0), p(1), minc(numext::real(m1.coeff(0))), maxc(numext::real(m1.coeff(0)));
  24. for(int j = 0; j < cols; j++)
  25. for(int i = 0; i < rows; i++)
  26. {
  27. s += m1(i,j);
  28. p *= m1_for_prod(i,j);
  29. minc = (std::min)(numext::real(minc), numext::real(m1(i,j)));
  30. maxc = (std::max)(numext::real(maxc), numext::real(m1(i,j)));
  31. }
  32. const Scalar mean = s/Scalar(RealScalar(rows*cols));
  33. VERIFY_IS_APPROX(m1.sum(), s);
  34. VERIFY_IS_APPROX(m1.mean(), mean);
  35. VERIFY_IS_APPROX(m1_for_prod.prod(), p);
  36. VERIFY_IS_APPROX(m1.real().minCoeff(), numext::real(minc));
  37. VERIFY_IS_APPROX(m1.real().maxCoeff(), numext::real(maxc));
  38. // test slice vectorization assuming assign is ok
  39. Index r0 = internal::random<Index>(0,rows-1);
  40. Index c0 = internal::random<Index>(0,cols-1);
  41. Index r1 = internal::random<Index>(r0+1,rows)-r0;
  42. Index c1 = internal::random<Index>(c0+1,cols)-c0;
  43. VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).sum(), m1.block(r0,c0,r1,c1).eval().sum());
  44. VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).mean(), m1.block(r0,c0,r1,c1).eval().mean());
  45. VERIFY_IS_APPROX(m1_for_prod.block(r0,c0,r1,c1).prod(), m1_for_prod.block(r0,c0,r1,c1).eval().prod());
  46. VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).real().minCoeff(), m1.block(r0,c0,r1,c1).real().eval().minCoeff());
  47. VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).real().maxCoeff(), m1.block(r0,c0,r1,c1).real().eval().maxCoeff());
  48. // test empty objects
  49. VERIFY_IS_APPROX(m1.block(r0,c0,0,0).sum(), Scalar(0));
  50. VERIFY_IS_APPROX(m1.block(r0,c0,0,0).prod(), Scalar(1));
  51. }
  52. template<typename VectorType> void vectorRedux(const VectorType& w)
  53. {
  54. using std::abs;
  55. typedef typename VectorType::Index Index;
  56. typedef typename VectorType::Scalar Scalar;
  57. typedef typename NumTraits<Scalar>::Real RealScalar;
  58. Index size = w.size();
  59. VectorType v = VectorType::Random(size);
  60. VectorType v_for_prod = VectorType::Ones(size) + Scalar(0.2) * v; // see comment above declaration of m1_for_prod
  61. for(int i = 1; i < size; i++)
  62. {
  63. Scalar s(0), p(1);
  64. RealScalar minc(numext::real(v.coeff(0))), maxc(numext::real(v.coeff(0)));
  65. for(int j = 0; j < i; j++)
  66. {
  67. s += v[j];
  68. p *= v_for_prod[j];
  69. minc = (std::min)(minc, numext::real(v[j]));
  70. maxc = (std::max)(maxc, numext::real(v[j]));
  71. }
  72. VERIFY_IS_MUCH_SMALLER_THAN(abs(s - v.head(i).sum()), Scalar(1));
  73. VERIFY_IS_APPROX(p, v_for_prod.head(i).prod());
  74. VERIFY_IS_APPROX(minc, v.real().head(i).minCoeff());
  75. VERIFY_IS_APPROX(maxc, v.real().head(i).maxCoeff());
  76. }
  77. for(int i = 0; i < size-1; i++)
  78. {
  79. Scalar s(0), p(1);
  80. RealScalar minc(numext::real(v.coeff(i))), maxc(numext::real(v.coeff(i)));
  81. for(int j = i; j < size; j++)
  82. {
  83. s += v[j];
  84. p *= v_for_prod[j];
  85. minc = (std::min)(minc, numext::real(v[j]));
  86. maxc = (std::max)(maxc, numext::real(v[j]));
  87. }
  88. VERIFY_IS_MUCH_SMALLER_THAN(abs(s - v.tail(size-i).sum()), Scalar(1));
  89. VERIFY_IS_APPROX(p, v_for_prod.tail(size-i).prod());
  90. VERIFY_IS_APPROX(minc, v.real().tail(size-i).minCoeff());
  91. VERIFY_IS_APPROX(maxc, v.real().tail(size-i).maxCoeff());
  92. }
  93. for(int i = 0; i < size/2; i++)
  94. {
  95. Scalar s(0), p(1);
  96. RealScalar minc(numext::real(v.coeff(i))), maxc(numext::real(v.coeff(i)));
  97. for(int j = i; j < size-i; j++)
  98. {
  99. s += v[j];
  100. p *= v_for_prod[j];
  101. minc = (std::min)(minc, numext::real(v[j]));
  102. maxc = (std::max)(maxc, numext::real(v[j]));
  103. }
  104. VERIFY_IS_MUCH_SMALLER_THAN(abs(s - v.segment(i, size-2*i).sum()), Scalar(1));
  105. VERIFY_IS_APPROX(p, v_for_prod.segment(i, size-2*i).prod());
  106. VERIFY_IS_APPROX(minc, v.real().segment(i, size-2*i).minCoeff());
  107. VERIFY_IS_APPROX(maxc, v.real().segment(i, size-2*i).maxCoeff());
  108. }
  109. // test empty objects
  110. VERIFY_IS_APPROX(v.head(0).sum(), Scalar(0));
  111. VERIFY_IS_APPROX(v.tail(0).prod(), Scalar(1));
  112. VERIFY_RAISES_ASSERT(v.head(0).mean());
  113. VERIFY_RAISES_ASSERT(v.head(0).minCoeff());
  114. VERIFY_RAISES_ASSERT(v.head(0).maxCoeff());
  115. }
  116. void test_redux()
  117. {
  118. // the max size cannot be too large, otherwise reduxion operations obviously generate large errors.
  119. int maxsize = (std::min)(100,EIGEN_TEST_MAX_SIZE);
  120. TEST_SET_BUT_UNUSED_VARIABLE(maxsize);
  121. for(int i = 0; i < g_repeat; i++) {
  122. CALL_SUBTEST_1( matrixRedux(Matrix<float, 1, 1>()) );
  123. CALL_SUBTEST_1( matrixRedux(Array<float, 1, 1>()) );
  124. CALL_SUBTEST_2( matrixRedux(Matrix2f()) );
  125. CALL_SUBTEST_2( matrixRedux(Array2f()) );
  126. CALL_SUBTEST_3( matrixRedux(Matrix4d()) );
  127. CALL_SUBTEST_3( matrixRedux(Array4d()) );
  128. CALL_SUBTEST_4( matrixRedux(MatrixXcf(internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  129. CALL_SUBTEST_4( matrixRedux(ArrayXXcf(internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  130. CALL_SUBTEST_5( matrixRedux(MatrixXd (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  131. CALL_SUBTEST_5( matrixRedux(ArrayXXd (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  132. CALL_SUBTEST_6( matrixRedux(MatrixXi (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  133. CALL_SUBTEST_6( matrixRedux(ArrayXXi (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
  134. }
  135. for(int i = 0; i < g_repeat; i++) {
  136. CALL_SUBTEST_7( vectorRedux(Vector4f()) );
  137. CALL_SUBTEST_7( vectorRedux(Array4f()) );
  138. CALL_SUBTEST_5( vectorRedux(VectorXd(internal::random<int>(1,maxsize))) );
  139. CALL_SUBTEST_5( vectorRedux(ArrayXd(internal::random<int>(1,maxsize))) );
  140. CALL_SUBTEST_8( vectorRedux(VectorXf(internal::random<int>(1,maxsize))) );
  141. CALL_SUBTEST_8( vectorRedux(ArrayXf(internal::random<int>(1,maxsize))) );
  142. }
  143. }