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