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