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 // this hack is needed to make this file compiles with -pedantic (gcc)
     12 #ifdef __GNUC__
     13 #define throw(X)
     14 #endif
     15 
     16 #ifdef __INTEL_COMPILER
     17   // disable "warning #76: argument to macro is empty" produced by the above hack
     18   #pragma warning disable 76
     19 #endif
     20 
     21 // discard stack allocation as that too bypasses malloc
     22 #define EIGEN_STACK_ALLOCATION_LIMIT 0
     23 // any heap allocation will raise an assert
     24 #define EIGEN_NO_MALLOC
     25 
     26 #include "main.h"
     27 #include <Eigen/Cholesky>
     28 #include <Eigen/Eigenvalues>
     29 #include <Eigen/LU>
     30 #include <Eigen/QR>
     31 #include <Eigen/SVD>
     32 
     33 template<typename MatrixType> void nomalloc(const MatrixType& m)
     34 {
     35   /* this test check no dynamic memory allocation are issued with fixed-size matrices
     36   */
     37   typedef typename MatrixType::Index Index;
     38   typedef typename MatrixType::Scalar Scalar;
     39 
     40   Index rows = m.rows();
     41   Index cols = m.cols();
     42 
     43   MatrixType m1 = MatrixType::Random(rows, cols),
     44              m2 = MatrixType::Random(rows, cols),
     45              m3(rows, cols);
     46 
     47   Scalar s1 = internal::random<Scalar>();
     48 
     49   Index r = internal::random<Index>(0, rows-1),
     50         c = internal::random<Index>(0, cols-1);
     51 
     52   VERIFY_IS_APPROX((m1+m2)*s1,              s1*m1+s1*m2);
     53   VERIFY_IS_APPROX((m1+m2)(r,c), (m1(r,c))+(m2(r,c)));
     54   VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix());
     55   VERIFY_IS_APPROX((m1*m1.transpose())*m2,  m1*(m1.transpose()*m2));
     56 
     57   m2.col(0).noalias() = m1 * m1.col(0);
     58   m2.col(0).noalias() -= m1.adjoint() * m1.col(0);
     59   m2.col(0).noalias() -= m1 * m1.row(0).adjoint();
     60   m2.col(0).noalias() -= m1.adjoint() * m1.row(0).adjoint();
     61 
     62   m2.row(0).noalias() = m1.row(0) * m1;
     63   m2.row(0).noalias() -= m1.row(0) * m1.adjoint();
     64   m2.row(0).noalias() -= m1.col(0).adjoint() * m1;
     65   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint();
     66   VERIFY_IS_APPROX(m2,m2);
     67 
     68   m2.col(0).noalias() = m1.template triangularView<Upper>() * m1.col(0);
     69   m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.col(0);
     70   m2.col(0).noalias() -= m1.template triangularView<Upper>() * m1.row(0).adjoint();
     71   m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.row(0).adjoint();
     72 
     73   m2.row(0).noalias() = m1.row(0) * m1.template triangularView<Upper>();
     74   m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template triangularView<Upper>();
     75   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template triangularView<Upper>();
     76   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template triangularView<Upper>();
     77   VERIFY_IS_APPROX(m2,m2);
     78 
     79   m2.col(0).noalias() = m1.template selfadjointView<Upper>() * m1.col(0);
     80   m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.col(0);
     81   m2.col(0).noalias() -= m1.template selfadjointView<Upper>() * m1.row(0).adjoint();
     82   m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.row(0).adjoint();
     83 
     84   m2.row(0).noalias() = m1.row(0) * m1.template selfadjointView<Upper>();
     85   m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template selfadjointView<Upper>();
     86   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template selfadjointView<Upper>();
     87   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template selfadjointView<Upper>();
     88   VERIFY_IS_APPROX(m2,m2);
     89 
     90   m2.template selfadjointView<Lower>().rankUpdate(m1.col(0),-1);
     91   m2.template selfadjointView<Lower>().rankUpdate(m1.row(0),-1);
     92 
     93   // The following fancy matrix-matrix products are not safe yet regarding static allocation
     94 //   m1 += m1.template triangularView<Upper>() * m2.col(;
     95 //   m1.template selfadjointView<Lower>().rankUpdate(m2);
     96 //   m1 += m1.template triangularView<Upper>() * m2;
     97 //   m1 += m1.template selfadjointView<Lower>() * m2;
     98 //   VERIFY_IS_APPROX(m1,m1);
     99 }
    100 
    101 template<typename Scalar>
    102 void ctms_decompositions()
    103 {
    104   const int maxSize = 16;
    105   const int size    = 12;
    106 
    107   typedef Eigen::Matrix<Scalar,
    108                         Eigen::Dynamic, Eigen::Dynamic,
    109                         0,
    110                         maxSize, maxSize> Matrix;
    111 
    112   typedef Eigen::Matrix<Scalar,
    113                         Eigen::Dynamic, 1,
    114                         0,
    115                         maxSize, 1> Vector;
    116 
    117   typedef Eigen::Matrix<std::complex<Scalar>,
    118                         Eigen::Dynamic, Eigen::Dynamic,
    119                         0,
    120                         maxSize, maxSize> ComplexMatrix;
    121 
    122   const Matrix A(Matrix::Random(size, size)), B(Matrix::Random(size, size));
    123   Matrix X(size,size);
    124   const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
    125   const Matrix saA = A.adjoint() * A;
    126   const Vector b(Vector::Random(size));
    127   Vector x(size);
    128 
    129   // Cholesky module
    130   Eigen::LLT<Matrix>  LLT;  LLT.compute(A);
    131   X = LLT.solve(B);
    132   x = LLT.solve(b);
    133   Eigen::LDLT<Matrix> LDLT; LDLT.compute(A);
    134   X = LDLT.solve(B);
    135   x = LDLT.solve(b);
    136 
    137   // Eigenvalues module
    138   Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp;        hessDecomp.compute(complexA);
    139   Eigen::ComplexSchur<ComplexMatrix>            cSchur(size);      cSchur.compute(complexA);
    140   Eigen::ComplexEigenSolver<ComplexMatrix>      cEigSolver;        cEigSolver.compute(complexA);
    141   Eigen::EigenSolver<Matrix>                    eigSolver;         eigSolver.compute(A);
    142   Eigen::SelfAdjointEigenSolver<Matrix>         saEigSolver(size); saEigSolver.compute(saA);
    143   Eigen::Tridiagonalization<Matrix>             tridiag;           tridiag.compute(saA);
    144 
    145   // LU module
    146   Eigen::PartialPivLU<Matrix> ppLU; ppLU.compute(A);
    147   X = ppLU.solve(B);
    148   x = ppLU.solve(b);
    149   Eigen::FullPivLU<Matrix>    fpLU; fpLU.compute(A);
    150   X = fpLU.solve(B);
    151   x = fpLU.solve(b);
    152 
    153   // QR module
    154   Eigen::HouseholderQR<Matrix>        hQR;  hQR.compute(A);
    155   X = hQR.solve(B);
    156   x = hQR.solve(b);
    157   Eigen::ColPivHouseholderQR<Matrix>  cpQR; cpQR.compute(A);
    158   X = cpQR.solve(B);
    159   x = cpQR.solve(b);
    160   Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A);
    161   // FIXME X = fpQR.solve(B);
    162   x = fpQR.solve(b);
    163 
    164   // SVD module
    165   Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A, ComputeFullU | ComputeFullV);
    166 }
    167 
    168 void test_nomalloc()
    169 {
    170   // check that our operator new is indeed called:
    171   VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3)));
    172   CALL_SUBTEST_1(nomalloc(Matrix<float, 1, 1>()) );
    173   CALL_SUBTEST_2(nomalloc(Matrix4d()) );
    174   CALL_SUBTEST_3(nomalloc(Matrix<float,32,32>()) );
    175 
    176   // Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
    177   CALL_SUBTEST_4(ctms_decompositions<float>());
    178 
    179 }
    180