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) 2006-2008 Benoit Jacob <jacob.benoit.1 (at) 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 #include <Eigen/QR>
     12 
     13 template<typename Derived1, typename Derived2>
     14 bool areNotApprox(const MatrixBase<Derived1>& m1, const MatrixBase<Derived2>& m2, typename Derived1::RealScalar epsilon = NumTraits<typename Derived1::RealScalar>::dummy_precision())
     15 {
     16   return !((m1-m2).cwiseAbs2().maxCoeff() < epsilon * epsilon
     17                           * (std::max)(m1.cwiseAbs2().maxCoeff(), m2.cwiseAbs2().maxCoeff()));
     18 }
     19 
     20 template<typename MatrixType> void product(const MatrixType& m)
     21 {
     22   /* this test covers the following files:
     23      Identity.h Product.h
     24   */
     25   typedef typename MatrixType::Index Index;
     26   typedef typename MatrixType::Scalar Scalar;
     27   typedef typename NumTraits<Scalar>::NonInteger NonInteger;
     28   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> RowVectorType;
     29   typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> ColVectorType;
     30   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> RowSquareMatrixType;
     31   typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> ColSquareMatrixType;
     32   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime,
     33                          MatrixType::Flags&RowMajorBit?ColMajor:RowMajor> OtherMajorMatrixType;
     34 
     35   Index rows = m.rows();
     36   Index cols = m.cols();
     37 
     38   // this test relies a lot on Random.h, and there's not much more that we can do
     39   // to test it, hence I consider that we will have tested Random.h
     40   MatrixType m1 = MatrixType::Random(rows, cols),
     41              m2 = MatrixType::Random(rows, cols),
     42              m3(rows, cols);
     43   RowSquareMatrixType
     44              identity = RowSquareMatrixType::Identity(rows, rows),
     45              square = RowSquareMatrixType::Random(rows, rows),
     46              res = RowSquareMatrixType::Random(rows, rows);
     47   ColSquareMatrixType
     48              square2 = ColSquareMatrixType::Random(cols, cols),
     49              res2 = ColSquareMatrixType::Random(cols, cols);
     50   RowVectorType v1 = RowVectorType::Random(rows);
     51   ColVectorType vc2 = ColVectorType::Random(cols), vcres(cols);
     52   OtherMajorMatrixType tm1 = m1;
     53 
     54   Scalar s1 = internal::random<Scalar>();
     55 
     56   Index r  = internal::random<Index>(0, rows-1),
     57         c  = internal::random<Index>(0, cols-1),
     58         c2 = internal::random<Index>(0, cols-1);
     59 
     60   // begin testing Product.h: only associativity for now
     61   // (we use Transpose.h but this doesn't count as a test for it)
     62   VERIFY_IS_APPROX((m1*m1.transpose())*m2,  m1*(m1.transpose()*m2));
     63   m3 = m1;
     64   m3 *= m1.transpose() * m2;
     65   VERIFY_IS_APPROX(m3,                      m1 * (m1.transpose()*m2));
     66   VERIFY_IS_APPROX(m3,                      m1 * (m1.transpose()*m2));
     67 
     68   // continue testing Product.h: distributivity
     69   VERIFY_IS_APPROX(square*(m1 + m2),        square*m1+square*m2);
     70   VERIFY_IS_APPROX(square*(m1 - m2),        square*m1-square*m2);
     71 
     72   // continue testing Product.h: compatibility with ScalarMultiple.h
     73   VERIFY_IS_APPROX(s1*(square*m1),          (s1*square)*m1);
     74   VERIFY_IS_APPROX(s1*(square*m1),          square*(m1*s1));
     75 
     76   // test Product.h together with Identity.h
     77   VERIFY_IS_APPROX(v1,                      identity*v1);
     78   VERIFY_IS_APPROX(v1.transpose(),          v1.transpose() * identity);
     79   // again, test operator() to check const-qualification
     80   VERIFY_IS_APPROX(MatrixType::Identity(rows, cols)(r,c), static_cast<Scalar>(r==c));
     81 
     82   if (rows!=cols)
     83      VERIFY_RAISES_ASSERT(m3 = m1*m1);
     84 
     85   // test the previous tests were not screwed up because operator* returns 0
     86   // (we use the more accurate default epsilon)
     87   if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
     88   {
     89     VERIFY(areNotApprox(m1.transpose()*m2,m2.transpose()*m1));
     90   }
     91 
     92   // test optimized operator+= path
     93   res = square;
     94   res.noalias() += m1 * m2.transpose();
     95   VERIFY_IS_APPROX(res, square + m1 * m2.transpose());
     96   if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
     97   {
     98     VERIFY(areNotApprox(res,square + m2 * m1.transpose()));
     99   }
    100   vcres = vc2;
    101   vcres.noalias() += m1.transpose() * v1;
    102   VERIFY_IS_APPROX(vcres, vc2 + m1.transpose() * v1);
    103 
    104   // test optimized operator-= path
    105   res = square;
    106   res.noalias() -= m1 * m2.transpose();
    107   VERIFY_IS_APPROX(res, square - (m1 * m2.transpose()));
    108   if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
    109   {
    110     VERIFY(areNotApprox(res,square - m2 * m1.transpose()));
    111   }
    112   vcres = vc2;
    113   vcres.noalias() -= m1.transpose() * v1;
    114   VERIFY_IS_APPROX(vcres, vc2 - m1.transpose() * v1);
    115 
    116   tm1 = m1;
    117   VERIFY_IS_APPROX(tm1.transpose() * v1, m1.transpose() * v1);
    118   VERIFY_IS_APPROX(v1.transpose() * tm1, v1.transpose() * m1);
    119 
    120   // test submatrix and matrix/vector product
    121   for (int i=0; i<rows; ++i)
    122     res.row(i) = m1.row(i) * m2.transpose();
    123   VERIFY_IS_APPROX(res, m1 * m2.transpose());
    124   // the other way round:
    125   for (int i=0; i<rows; ++i)
    126     res.col(i) = m1 * m2.transpose().col(i);
    127   VERIFY_IS_APPROX(res, m1 * m2.transpose());
    128 
    129   res2 = square2;
    130   res2.noalias() += m1.transpose() * m2;
    131   VERIFY_IS_APPROX(res2, square2 + m1.transpose() * m2);
    132   if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
    133   {
    134     VERIFY(areNotApprox(res2,square2 + m2.transpose() * m1));
    135   }
    136 
    137   VERIFY_IS_APPROX(res.col(r).noalias() = square.adjoint() * square.col(r), (square.adjoint() * square.col(r)).eval());
    138   VERIFY_IS_APPROX(res.col(r).noalias() = square * square.col(r), (square * square.col(r)).eval());
    139 
    140   // inner product
    141   Scalar x = square2.row(c) * square2.col(c2);
    142   VERIFY_IS_APPROX(x, square2.row(c).transpose().cwiseProduct(square2.col(c2)).sum());
    143 }
    144