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 // discard stack allocation as that too bypasses malloc
12 #define EIGEN_STACK_ALLOCATION_LIMIT 0
13 // heap allocation will raise an assert if enabled at runtime
14 #define EIGEN_RUNTIME_NO_MALLOC
15 
16 #include "main.h"
17 #include <Eigen/Cholesky>
18 #include <Eigen/Eigenvalues>
19 #include <Eigen/LU>
20 #include <Eigen/QR>
21 #include <Eigen/SVD>
22 
nomalloc(const MatrixType & m)23 template<typename MatrixType> void nomalloc(const MatrixType& m)
24 {
25   /* this test check no dynamic memory allocation are issued with fixed-size matrices
26   */
27   typedef typename MatrixType::Index Index;
28   typedef typename MatrixType::Scalar Scalar;
29 
30   Index rows = m.rows();
31   Index cols = m.cols();
32 
33   MatrixType m1 = MatrixType::Random(rows, cols),
34              m2 = MatrixType::Random(rows, cols),
35              m3(rows, cols);
36 
37   Scalar s1 = internal::random<Scalar>();
38 
39   Index r = internal::random<Index>(0, rows-1),
40         c = internal::random<Index>(0, cols-1);
41 
42   VERIFY_IS_APPROX((m1+m2)*s1,              s1*m1+s1*m2);
43   VERIFY_IS_APPROX((m1+m2)(r,c), (m1(r,c))+(m2(r,c)));
44   VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix());
45   VERIFY_IS_APPROX((m1*m1.transpose())*m2,  m1*(m1.transpose()*m2));
46 
47   m2.col(0).noalias() = m1 * m1.col(0);
48   m2.col(0).noalias() -= m1.adjoint() * m1.col(0);
49   m2.col(0).noalias() -= m1 * m1.row(0).adjoint();
50   m2.col(0).noalias() -= m1.adjoint() * m1.row(0).adjoint();
51 
52   m2.row(0).noalias() = m1.row(0) * m1;
53   m2.row(0).noalias() -= m1.row(0) * m1.adjoint();
54   m2.row(0).noalias() -= m1.col(0).adjoint() * m1;
55   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint();
56   VERIFY_IS_APPROX(m2,m2);
57 
58   m2.col(0).noalias() = m1.template triangularView<Upper>() * m1.col(0);
59   m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.col(0);
60   m2.col(0).noalias() -= m1.template triangularView<Upper>() * m1.row(0).adjoint();
61   m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.row(0).adjoint();
62 
63   m2.row(0).noalias() = m1.row(0) * m1.template triangularView<Upper>();
64   m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template triangularView<Upper>();
65   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template triangularView<Upper>();
66   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template triangularView<Upper>();
67   VERIFY_IS_APPROX(m2,m2);
68 
69   m2.col(0).noalias() = m1.template selfadjointView<Upper>() * m1.col(0);
70   m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.col(0);
71   m2.col(0).noalias() -= m1.template selfadjointView<Upper>() * m1.row(0).adjoint();
72   m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.row(0).adjoint();
73 
74   m2.row(0).noalias() = m1.row(0) * m1.template selfadjointView<Upper>();
75   m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template selfadjointView<Upper>();
76   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template selfadjointView<Upper>();
77   m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template selfadjointView<Upper>();
78   VERIFY_IS_APPROX(m2,m2);
79 
80   m2.template selfadjointView<Lower>().rankUpdate(m1.col(0),-1);
81   m2.template selfadjointView<Upper>().rankUpdate(m1.row(0),-1);
82   m2.template selfadjointView<Lower>().rankUpdate(m1.col(0), m1.col(0)); // rank-2
83 
84   // The following fancy matrix-matrix products are not safe yet regarding static allocation
85   m2.template selfadjointView<Lower>().rankUpdate(m1);
86   m2 += m2.template triangularView<Upper>() * m1;
87   m2.template triangularView<Upper>() = m2 * m2;
88   m1 += m1.template selfadjointView<Lower>() * m2;
89   VERIFY_IS_APPROX(m2,m2);
90 }
91 
92 template<typename Scalar>
ctms_decompositions()93 void ctms_decompositions()
94 {
95   const int maxSize = 16;
96   const int size    = 12;
97 
98   typedef Eigen::Matrix<Scalar,
99                         Eigen::Dynamic, Eigen::Dynamic,
100                         0,
101                         maxSize, maxSize> Matrix;
102 
103   typedef Eigen::Matrix<Scalar,
104                         Eigen::Dynamic, 1,
105                         0,
106                         maxSize, 1> Vector;
107 
108   typedef Eigen::Matrix<std::complex<Scalar>,
109                         Eigen::Dynamic, Eigen::Dynamic,
110                         0,
111                         maxSize, maxSize> ComplexMatrix;
112 
113   const Matrix A(Matrix::Random(size, size)), B(Matrix::Random(size, size));
114   Matrix X(size,size);
115   const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
116   const Matrix saA = A.adjoint() * A;
117   const Vector b(Vector::Random(size));
118   Vector x(size);
119 
120   // Cholesky module
121   Eigen::LLT<Matrix>  LLT;  LLT.compute(A);
122   X = LLT.solve(B);
123   x = LLT.solve(b);
124   Eigen::LDLT<Matrix> LDLT; LDLT.compute(A);
125   X = LDLT.solve(B);
126   x = LDLT.solve(b);
127 
128   // Eigenvalues module
129   Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp;        hessDecomp.compute(complexA);
130   Eigen::ComplexSchur<ComplexMatrix>            cSchur(size);      cSchur.compute(complexA);
131   Eigen::ComplexEigenSolver<ComplexMatrix>      cEigSolver;        cEigSolver.compute(complexA);
132   Eigen::EigenSolver<Matrix>                    eigSolver;         eigSolver.compute(A);
133   Eigen::SelfAdjointEigenSolver<Matrix>         saEigSolver(size); saEigSolver.compute(saA);
134   Eigen::Tridiagonalization<Matrix>             tridiag;           tridiag.compute(saA);
135 
136   // LU module
137   Eigen::PartialPivLU<Matrix> ppLU; ppLU.compute(A);
138   X = ppLU.solve(B);
139   x = ppLU.solve(b);
140   Eigen::FullPivLU<Matrix>    fpLU; fpLU.compute(A);
141   X = fpLU.solve(B);
142   x = fpLU.solve(b);
143 
144   // QR module
145   Eigen::HouseholderQR<Matrix>        hQR;  hQR.compute(A);
146   X = hQR.solve(B);
147   x = hQR.solve(b);
148   Eigen::ColPivHouseholderQR<Matrix>  cpQR; cpQR.compute(A);
149   X = cpQR.solve(B);
150   x = cpQR.solve(b);
151   Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A);
152   // FIXME X = fpQR.solve(B);
153   x = fpQR.solve(b);
154 
155   // SVD module
156   Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A, ComputeFullU | ComputeFullV);
157 }
158 
test_zerosized()159 void test_zerosized() {
160   // default constructors:
161   Eigen::MatrixXd A;
162   Eigen::VectorXd v;
163   // explicit zero-sized:
164   Eigen::ArrayXXd A0(0,0);
165   Eigen::ArrayXd v0(0);
166 
167   // assigning empty objects to each other:
168   A=A0;
169   v=v0;
170 }
171 
test_reference(const MatrixType & m)172 template<typename MatrixType> void test_reference(const MatrixType& m) {
173   typedef typename MatrixType::Scalar Scalar;
174   enum { Flag          =  MatrixType::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor};
175   enum { TransposeFlag = !MatrixType::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor};
176   typename MatrixType::Index rows = m.rows(), cols=m.cols();
177   typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flag         > MatrixX;
178   typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, TransposeFlag> MatrixXT;
179   // Dynamic reference:
180   typedef Eigen::Ref<const MatrixX  > Ref;
181   typedef Eigen::Ref<const MatrixXT > RefT;
182 
183   Ref r1(m);
184   Ref r2(m.block(rows/3, cols/4, rows/2, cols/2));
185   RefT r3(m.transpose());
186   RefT r4(m.topLeftCorner(rows/2, cols/2).transpose());
187 
188   VERIFY_RAISES_ASSERT(RefT r5(m));
189   VERIFY_RAISES_ASSERT(Ref r6(m.transpose()));
190   VERIFY_RAISES_ASSERT(Ref r7(Scalar(2) * m));
191 
192   // Copy constructors shall also never malloc
193   Ref r8 = r1;
194   RefT r9 = r3;
195 
196   // Initializing from a compatible Ref shall also never malloc
197   Eigen::Ref<const MatrixX, Unaligned, Stride<Dynamic, Dynamic> > r10=r8, r11=m;
198 
199   // Initializing from an incompatible Ref will malloc:
200   typedef Eigen::Ref<const MatrixX, Aligned> RefAligned;
201   VERIFY_RAISES_ASSERT(RefAligned r12=r10);
202   VERIFY_RAISES_ASSERT(Ref r13=r10); // r10 has more dynamic strides
203 
204 }
205 
test_nomalloc()206 void test_nomalloc()
207 {
208   // create some dynamic objects
209   Eigen::MatrixXd M1 = MatrixXd::Random(3,3);
210   Ref<const MatrixXd> R1 = 2.0*M1; // Ref requires temporary
211 
212   // from here on prohibit malloc:
213   Eigen::internal::set_is_malloc_allowed(false);
214 
215   // check that our operator new is indeed called:
216   VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3)));
217   CALL_SUBTEST_1(nomalloc(Matrix<float, 1, 1>()) );
218   CALL_SUBTEST_2(nomalloc(Matrix4d()) );
219   CALL_SUBTEST_3(nomalloc(Matrix<float,32,32>()) );
220 
221   // Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
222   CALL_SUBTEST_4(ctms_decompositions<float>());
223 
224   CALL_SUBTEST_5(test_zerosized());
225 
226   CALL_SUBTEST_6(test_reference(Matrix<float,32,32>()));
227   CALL_SUBTEST_7(test_reference(R1));
228   CALL_SUBTEST_8(Ref<MatrixXd> R2 = M1.topRows<2>(); test_reference(R2));
229 }
230