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 Gael Guennebaud <gael.guennebaud (at) inria.fr>
      5 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1 (at) gmail.com>
      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 EIGEN2_SUPPORT
     12 #define EIGEN_NO_STATIC_ASSERT
     13 #include "main.h"
     14 #include <functional>
     15 
     16 #ifdef min
     17 #undef min
     18 #endif
     19 
     20 #ifdef max
     21 #undef max
     22 #endif
     23 
     24 using namespace std;
     25 
     26 template<typename Scalar> struct AddIfNull {
     27     const Scalar operator() (const Scalar a, const Scalar b) const {return a<=1e-3 ? b : a;}
     28     enum { Cost = NumTraits<Scalar>::AddCost };
     29 };
     30 
     31 template<typename MatrixType>
     32 typename Eigen::internal::enable_if<!NumTraits<typename MatrixType::Scalar>::IsInteger,typename MatrixType::Scalar>::type
     33 cwiseops_real_only(MatrixType& m1, MatrixType& m2, MatrixType& m3, MatrixType& mones)
     34 {
     35   typedef typename MatrixType::Scalar Scalar;
     36   typedef typename NumTraits<Scalar>::Real RealScalar;
     37 
     38   VERIFY_IS_APPROX(m1.cwise() / m2,    m1.cwise() * (m2.cwise().inverse()));
     39   m3 = m1.cwise().abs().cwise().sqrt();
     40   VERIFY_IS_APPROX(m3.cwise().square(), m1.cwise().abs());
     41   VERIFY_IS_APPROX(m1.cwise().square().cwise().sqrt(), m1.cwise().abs());
     42   VERIFY_IS_APPROX(m1.cwise().abs().cwise().log().cwise().exp() , m1.cwise().abs());
     43 
     44   VERIFY_IS_APPROX(m1.cwise().pow(2), m1.cwise().square());
     45   m3 = (m1.cwise().abs().cwise()<=RealScalar(0.01)).select(mones,m1);
     46   VERIFY_IS_APPROX(m3.cwise().pow(-1), m3.cwise().inverse());
     47   m3 = m1.cwise().abs();
     48   VERIFY_IS_APPROX(m3.cwise().pow(RealScalar(0.5)), m3.cwise().sqrt());
     49 
     50 //   VERIFY_IS_APPROX(m1.cwise().tan(), m1.cwise().sin().cwise() / m1.cwise().cos());
     51   VERIFY_IS_APPROX(mones, m1.cwise().sin().cwise().square() + m1.cwise().cos().cwise().square());
     52   m3 = m1;
     53   m3.cwise() /= m2;
     54   VERIFY_IS_APPROX(m3, m1.cwise() / m2);
     55 
     56   return Scalar(0);
     57 }
     58 
     59 template<typename MatrixType>
     60 typename Eigen::internal::enable_if<NumTraits<typename MatrixType::Scalar>::IsInteger,typename MatrixType::Scalar>::type
     61 cwiseops_real_only(MatrixType& , MatrixType& , MatrixType& , MatrixType& )
     62 {
     63   return 0;
     64 }
     65 
     66 template<typename MatrixType> void cwiseops(const MatrixType& m)
     67 {
     68   typedef typename MatrixType::Index Index;
     69   typedef typename MatrixType::Scalar Scalar;
     70   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
     71 
     72   Index rows = m.rows();
     73   Index cols = m.cols();
     74 
     75   MatrixType m1 = MatrixType::Random(rows, cols),
     76              m1bis = m1,
     77              m2 = MatrixType::Random(rows, cols),
     78              m3(rows, cols),
     79              m4(rows, cols),
     80              mzero = MatrixType::Zero(rows, cols),
     81              mones = MatrixType::Ones(rows, cols),
     82              identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
     83                               ::Identity(rows, rows);
     84   VectorType vzero = VectorType::Zero(rows),
     85              vones = VectorType::Ones(rows),
     86              v3(rows);
     87 
     88   Index r = internal::random<Index>(0, rows-1),
     89         c = internal::random<Index>(0, cols-1);
     90 
     91   Scalar s1 = internal::random<Scalar>();
     92 
     93   // test Zero, Ones, Constant, and the set* variants
     94   m3 = MatrixType::Constant(rows, cols, s1);
     95   for (int j=0; j<cols; ++j)
     96     for (int i=0; i<rows; ++i)
     97     {
     98       VERIFY_IS_APPROX(mzero(i,j), Scalar(0));
     99       VERIFY_IS_APPROX(mones(i,j), Scalar(1));
    100       VERIFY_IS_APPROX(m3(i,j), s1);
    101     }
    102   VERIFY(mzero.isZero());
    103   VERIFY(mones.isOnes());
    104   VERIFY(m3.isConstant(s1));
    105   VERIFY(identity.isIdentity());
    106   VERIFY_IS_APPROX(m4.setConstant(s1), m3);
    107   VERIFY_IS_APPROX(m4.setConstant(rows,cols,s1), m3);
    108   VERIFY_IS_APPROX(m4.setZero(), mzero);
    109   VERIFY_IS_APPROX(m4.setZero(rows,cols), mzero);
    110   VERIFY_IS_APPROX(m4.setOnes(), mones);
    111   VERIFY_IS_APPROX(m4.setOnes(rows,cols), mones);
    112   m4.fill(s1);
    113   VERIFY_IS_APPROX(m4, m3);
    114 
    115   VERIFY_IS_APPROX(v3.setConstant(rows, s1), VectorType::Constant(rows,s1));
    116   VERIFY_IS_APPROX(v3.setZero(rows), vzero);
    117   VERIFY_IS_APPROX(v3.setOnes(rows), vones);
    118 
    119   m2 = m2.template binaryExpr<AddIfNull<Scalar> >(mones);
    120 
    121   VERIFY_IS_APPROX(m1.cwise().pow(2), m1.cwise().abs2());
    122   VERIFY_IS_APPROX(m1.cwise().pow(2), m1.cwise().square());
    123   VERIFY_IS_APPROX(m1.cwise().pow(3), m1.cwise().cube());
    124 
    125   VERIFY_IS_APPROX(m1 + mones, m1.cwise()+Scalar(1));
    126   VERIFY_IS_APPROX(m1 - mones, m1.cwise()-Scalar(1));
    127   m3 = m1; m3.cwise() += 1;
    128   VERIFY_IS_APPROX(m1 + mones, m3);
    129   m3 = m1; m3.cwise() -= 1;
    130   VERIFY_IS_APPROX(m1 - mones, m3);
    131 
    132   VERIFY_IS_APPROX(m2, m2.cwise() * mones);
    133   VERIFY_IS_APPROX(m1.cwise() * m2,  m2.cwise() * m1);
    134   m3 = m1;
    135   m3.cwise() *= m2;
    136   VERIFY_IS_APPROX(m3, m1.cwise() * m2);
    137 
    138   VERIFY_IS_APPROX(mones,    m2.cwise()/m2);
    139 
    140   // check min
    141   VERIFY_IS_APPROX( m1.cwise().min(m2), m2.cwise().min(m1) );
    142   VERIFY_IS_APPROX( m1.cwise().min(m1+mones), m1 );
    143   VERIFY_IS_APPROX( m1.cwise().min(m1-mones), m1-mones );
    144 
    145   // check max
    146   VERIFY_IS_APPROX( m1.cwise().max(m2), m2.cwise().max(m1) );
    147   VERIFY_IS_APPROX( m1.cwise().max(m1-mones), m1 );
    148   VERIFY_IS_APPROX( m1.cwise().max(m1+mones), m1+mones );
    149 
    150   VERIFY( (m1.cwise() == m1).all() );
    151   VERIFY( (m1.cwise() != m2).any() );
    152   VERIFY(!(m1.cwise() == (m1+mones)).any() );
    153   if (rows*cols>1)
    154   {
    155     m3 = m1;
    156     m3(r,c) += 1;
    157     VERIFY( (m1.cwise() == m3).any() );
    158     VERIFY( !(m1.cwise() == m3).all() );
    159   }
    160   VERIFY( (m1.cwise().min(m2).cwise() <= m2).all() );
    161   VERIFY( (m1.cwise().max(m2).cwise() >= m2).all() );
    162   VERIFY( (m1.cwise().min(m2).cwise() < (m1+mones)).all() );
    163   VERIFY( (m1.cwise().max(m2).cwise() > (m1-mones)).all() );
    164 
    165   VERIFY( (m1.cwise()<m1.unaryExpr(bind2nd(plus<Scalar>(), Scalar(1)))).all() );
    166   VERIFY( !(m1.cwise()<m1bis.unaryExpr(bind2nd(minus<Scalar>(), Scalar(1)))).all() );
    167   VERIFY( !(m1.cwise()>m1bis.unaryExpr(bind2nd(plus<Scalar>(), Scalar(1)))).any() );
    168 
    169   cwiseops_real_only(m1, m2, m3, mones);
    170 }
    171 
    172 void test_cwiseop()
    173 {
    174   for(int i = 0; i < g_repeat ; i++) {
    175     CALL_SUBTEST_1( cwiseops(Matrix<float, 1, 1>()) );
    176     CALL_SUBTEST_2( cwiseops(Matrix4d()) );
    177     CALL_SUBTEST_3( cwiseops(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
    178     CALL_SUBTEST_4( cwiseops(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
    179     CALL_SUBTEST_5( cwiseops(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
    180     CALL_SUBTEST_6( cwiseops(MatrixXd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
    181   }
    182 }
    183