Home | History | Annotate | Download | only in test
      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 (at) gmail.com>
      5 // Copyright (C) 2015 Gael Guennebaud <gael.guennebaud (at) 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 
     11 #define TEST_ENABLE_TEMPORARY_TRACKING
     12 
     13 #include "main.h"
     14 
     15 template<typename MatrixType> void matrixRedux(const MatrixType& m)
     16 {
     17   typedef typename MatrixType::Index Index;
     18   typedef typename MatrixType::Scalar Scalar;
     19   typedef typename MatrixType::RealScalar RealScalar;
     20 
     21   Index rows = m.rows();
     22   Index cols = m.cols();
     23 
     24   MatrixType m1 = MatrixType::Random(rows, cols);
     25 
     26   // The entries of m1 are uniformly distributed in [0,1], so m1.prod() is very small. This may lead to test
     27   // failures if we underflow into denormals. Thus, we scale so that entries are close to 1.
     28   MatrixType m1_for_prod = MatrixType::Ones(rows, cols) + RealScalar(0.2) * m1;
     29 
     30   VERIFY_IS_MUCH_SMALLER_THAN(MatrixType::Zero(rows, cols).sum(), Scalar(1));
     31   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
     32   Scalar s(0), p(1), minc(numext::real(m1.coeff(0))), maxc(numext::real(m1.coeff(0)));
     33   for(int j = 0; j < cols; j++)
     34   for(int i = 0; i < rows; i++)
     35   {
     36     s += m1(i,j);
     37     p *= m1_for_prod(i,j);
     38     minc = (std::min)(numext::real(minc), numext::real(m1(i,j)));
     39     maxc = (std::max)(numext::real(maxc), numext::real(m1(i,j)));
     40   }
     41   const Scalar mean = s/Scalar(RealScalar(rows*cols));
     42 
     43   VERIFY_IS_APPROX(m1.sum(), s);
     44   VERIFY_IS_APPROX(m1.mean(), mean);
     45   VERIFY_IS_APPROX(m1_for_prod.prod(), p);
     46   VERIFY_IS_APPROX(m1.real().minCoeff(), numext::real(minc));
     47   VERIFY_IS_APPROX(m1.real().maxCoeff(), numext::real(maxc));
     48 
     49   // test slice vectorization assuming assign is ok
     50   Index r0 = internal::random<Index>(0,rows-1);
     51   Index c0 = internal::random<Index>(0,cols-1);
     52   Index r1 = internal::random<Index>(r0+1,rows)-r0;
     53   Index c1 = internal::random<Index>(c0+1,cols)-c0;
     54   VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).sum(), m1.block(r0,c0,r1,c1).eval().sum());
     55   VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).mean(), m1.block(r0,c0,r1,c1).eval().mean());
     56   VERIFY_IS_APPROX(m1_for_prod.block(r0,c0,r1,c1).prod(), m1_for_prod.block(r0,c0,r1,c1).eval().prod());
     57   VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).real().minCoeff(), m1.block(r0,c0,r1,c1).real().eval().minCoeff());
     58   VERIFY_IS_APPROX(m1.block(r0,c0,r1,c1).real().maxCoeff(), m1.block(r0,c0,r1,c1).real().eval().maxCoeff());
     59 
     60   // regression for bug 1090
     61   const int R1 = MatrixType::RowsAtCompileTime>=2 ? MatrixType::RowsAtCompileTime/2 : 6;
     62   const int C1 = MatrixType::ColsAtCompileTime>=2 ? MatrixType::ColsAtCompileTime/2 : 6;
     63   if(R1<=rows-r0 && C1<=cols-c0)
     64   {
     65     VERIFY_IS_APPROX( (m1.template block<R1,C1>(r0,c0).sum()), m1.block(r0,c0,R1,C1).sum() );
     66   }
     67 
     68   // test empty objects
     69   VERIFY_IS_APPROX(m1.block(r0,c0,0,0).sum(),   Scalar(0));
     70   VERIFY_IS_APPROX(m1.block(r0,c0,0,0).prod(),  Scalar(1));
     71 
     72   // test nesting complex expression
     73   VERIFY_EVALUATION_COUNT( (m1.matrix()*m1.matrix().transpose()).sum(), (MatrixType::IsVectorAtCompileTime && MatrixType::SizeAtCompileTime!=1 ? 0 : 1) );
     74   Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> m2(rows,rows);
     75   m2.setRandom();
     76   VERIFY_EVALUATION_COUNT( ((m1.matrix()*m1.matrix().transpose())+m2).sum(),(MatrixType::IsVectorAtCompileTime && MatrixType::SizeAtCompileTime!=1 ? 0 : 1));
     77 }
     78 
     79 template<typename VectorType> void vectorRedux(const VectorType& w)
     80 {
     81   using std::abs;
     82   typedef typename VectorType::Index Index;
     83   typedef typename VectorType::Scalar Scalar;
     84   typedef typename NumTraits<Scalar>::Real RealScalar;
     85   Index size = w.size();
     86 
     87   VectorType v = VectorType::Random(size);
     88   VectorType v_for_prod = VectorType::Ones(size) + Scalar(0.2) * v; // see comment above declaration of m1_for_prod
     89 
     90   for(int i = 1; i < size; i++)
     91   {
     92     Scalar s(0), p(1);
     93     RealScalar minc(numext::real(v.coeff(0))), maxc(numext::real(v.coeff(0)));
     94     for(int j = 0; j < i; 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.head(i).sum()), Scalar(1));
    102     VERIFY_IS_APPROX(p, v_for_prod.head(i).prod());
    103     VERIFY_IS_APPROX(minc, v.real().head(i).minCoeff());
    104     VERIFY_IS_APPROX(maxc, v.real().head(i).maxCoeff());
    105   }
    106 
    107   for(int i = 0; i < size-1; 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; 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.tail(size-i).sum()), Scalar(1));
    119     VERIFY_IS_APPROX(p, v_for_prod.tail(size-i).prod());
    120     VERIFY_IS_APPROX(minc, v.real().tail(size-i).minCoeff());
    121     VERIFY_IS_APPROX(maxc, v.real().tail(size-i).maxCoeff());
    122   }
    123 
    124   for(int i = 0; i < size/2; i++)
    125   {
    126     Scalar s(0), p(1);
    127     RealScalar minc(numext::real(v.coeff(i))), maxc(numext::real(v.coeff(i)));
    128     for(int j = i; j < size-i; j++)
    129     {
    130       s += v[j];
    131       p *= v_for_prod[j];
    132       minc = (std::min)(minc, numext::real(v[j]));
    133       maxc = (std::max)(maxc, numext::real(v[j]));
    134     }
    135     VERIFY_IS_MUCH_SMALLER_THAN(abs(s - v.segment(i, size-2*i).sum()), Scalar(1));
    136     VERIFY_IS_APPROX(p, v_for_prod.segment(i, size-2*i).prod());
    137     VERIFY_IS_APPROX(minc, v.real().segment(i, size-2*i).minCoeff());
    138     VERIFY_IS_APPROX(maxc, v.real().segment(i, size-2*i).maxCoeff());
    139   }
    140 
    141   // test empty objects
    142   VERIFY_IS_APPROX(v.head(0).sum(),   Scalar(0));
    143   VERIFY_IS_APPROX(v.tail(0).prod(),  Scalar(1));
    144   VERIFY_RAISES_ASSERT(v.head(0).mean());
    145   VERIFY_RAISES_ASSERT(v.head(0).minCoeff());
    146   VERIFY_RAISES_ASSERT(v.head(0).maxCoeff());
    147 }
    148 
    149 void test_redux()
    150 {
    151   // the max size cannot be too large, otherwise reduxion operations obviously generate large errors.
    152   int maxsize = (std::min)(100,EIGEN_TEST_MAX_SIZE);
    153   TEST_SET_BUT_UNUSED_VARIABLE(maxsize);
    154   for(int i = 0; i < g_repeat; i++) {
    155     CALL_SUBTEST_1( matrixRedux(Matrix<float, 1, 1>()) );
    156     CALL_SUBTEST_1( matrixRedux(Array<float, 1, 1>()) );
    157     CALL_SUBTEST_2( matrixRedux(Matrix2f()) );
    158     CALL_SUBTEST_2( matrixRedux(Array2f()) );
    159     CALL_SUBTEST_2( matrixRedux(Array22f()) );
    160     CALL_SUBTEST_3( matrixRedux(Matrix4d()) );
    161     CALL_SUBTEST_3( matrixRedux(Array4d()) );
    162     CALL_SUBTEST_3( matrixRedux(Array44d()) );
    163     CALL_SUBTEST_4( matrixRedux(MatrixXcf(internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
    164     CALL_SUBTEST_4( matrixRedux(ArrayXXcf(internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
    165     CALL_SUBTEST_5( matrixRedux(MatrixXd (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
    166     CALL_SUBTEST_5( matrixRedux(ArrayXXd (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
    167     CALL_SUBTEST_6( matrixRedux(MatrixXi (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
    168     CALL_SUBTEST_6( matrixRedux(ArrayXXi (internal::random<int>(1,maxsize), internal::random<int>(1,maxsize))) );
    169   }
    170   for(int i = 0; i < g_repeat; i++) {
    171     CALL_SUBTEST_7( vectorRedux(Vector4f()) );
    172     CALL_SUBTEST_7( vectorRedux(Array4f()) );
    173     CALL_SUBTEST_5( vectorRedux(VectorXd(internal::random<int>(1,maxsize))) );
    174     CALL_SUBTEST_5( vectorRedux(ArrayXd(internal::random<int>(1,maxsize))) );
    175     CALL_SUBTEST_8( vectorRedux(VectorXf(internal::random<int>(1,maxsize))) );
    176     CALL_SUBTEST_8( vectorRedux(ArrayXf(internal::random<int>(1,maxsize))) );
    177   }
    178 }
    179