blob: 30666421099b67ec3a596d55b299906cd26de2a1 [file] [log] [blame]
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +00003//
Gael Guennebaud28e64b02010-06-24 23:21:58 +02004// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
Benoit Jacob00f89a82008-11-24 13:40:43 +00005// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +00006//
Benoit Jacob69124cf2012-07-13 14:42:47 -04007// 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/.
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000010
Gael Guennebaude3fac692008-06-08 15:03:23 +000011// this hack is needed to make this file compiles with -pedantic (gcc)
Benoit Jacob9e3c7312009-01-22 00:09:34 +000012#ifdef __GNUC__
Gael Guennebaud60726f92008-06-04 10:15:48 +000013#define throw(X)
Benoit Jacob9e3c7312009-01-22 00:09:34 +000014#endif
Gael Guennebaud61a29952013-02-26 18:10:19 +010015
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
Benoit Jacob1d52bd42009-01-08 15:20:21 +000021// discard stack allocation as that too bypasses malloc
Benoit Jacobd6712052008-12-31 00:25:31 +000022#define EIGEN_STACK_ALLOCATION_LIMIT 0
Benoit Jacob1d52bd42009-01-08 15:20:21 +000023// any heap allocation will raise an assert
Christoph Hertzberg12d59462014-09-30 14:57:54 +020024#define EIGEN_RUNTIME_NO_MALLOC
Benoit Jacobd6712052008-12-31 00:25:31 +000025
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000026#include "main.h"
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +010027#include <Eigen/Cholesky>
28#include <Eigen/Eigenvalues>
29#include <Eigen/LU>
30#include <Eigen/QR>
31#include <Eigen/SVD>
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000032
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000033template<typename MatrixType> void nomalloc(const MatrixType& m)
34{
35 /* this test check no dynamic memory allocation are issued with fixed-size matrices
36 */
Hauke Heibelf1679c72010-06-20 17:37:56 +020037 typedef typename MatrixType::Index Index;
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000038 typedef typename MatrixType::Scalar Scalar;
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000039
Hauke Heibelf1679c72010-06-20 17:37:56 +020040 Index rows = m.rows();
41 Index cols = m.cols();
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000042
Gael Guennebaudc10f0692008-07-21 00:34:46 +000043 MatrixType m1 = MatrixType::Random(rows, cols),
44 m2 = MatrixType::Random(rows, cols),
Benoit Jacob0609dbe2011-10-31 00:51:36 -040045 m3(rows, cols);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000046
Benoit Jacob47160402010-10-25 10:15:22 -040047 Scalar s1 = internal::random<Scalar>();
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000048
Benoit Jacob47160402010-10-25 10:15:22 -040049 Index r = internal::random<Index>(0, rows-1),
50 c = internal::random<Index>(0, cols-1);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000051
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)));
Gael Guennebaud71a171c2010-01-04 19:13:08 +010054 VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix());
Gael Guennebaud3d8e1792011-01-27 18:02:49 +010055 VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
Gael Guennebaudc60818f2011-02-01 11:38:46 +010056
Gael Guennebaudf46ace62011-01-31 21:30:27 +010057 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();
Gael Guennebaudc60818f2011-02-01 11:38:46 +010066 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);
Gael Guennebaud59af20b2011-02-01 16:46:35 +010089
90 m2.template selfadjointView<Lower>().rankUpdate(m1.col(0),-1);
91 m2.template selfadjointView<Lower>().rankUpdate(m1.row(0),-1);
Gael Guennebaudc60818f2011-02-01 11:38:46 +010092
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);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000099}
100
Gael Guennebaud8a380472010-07-04 15:35:21 +0200101template<typename Scalar>
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100102void ctms_decompositions()
103{
104 const int maxSize = 16;
105 const int size = 12;
106
Gael Guennebaud8a380472010-07-04 15:35:21 +0200107 typedef Eigen::Matrix<Scalar,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100108 Eigen::Dynamic, Eigen::Dynamic,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500109 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100110 maxSize, maxSize> Matrix;
111
Gael Guennebaud8a380472010-07-04 15:35:21 +0200112 typedef Eigen::Matrix<Scalar,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100113 Eigen::Dynamic, 1,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500114 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100115 maxSize, 1> Vector;
116
Gael Guennebaud8a380472010-07-04 15:35:21 +0200117 typedef Eigen::Matrix<std::complex<Scalar>,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100118 Eigen::Dynamic, Eigen::Dynamic,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500119 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100120 maxSize, maxSize> ComplexMatrix;
121
Gael Guennebaud88e05102012-06-12 13:12:47 +0200122 const Matrix A(Matrix::Random(size, size)), B(Matrix::Random(size, size));
123 Matrix X(size,size);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100124 const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
Gael Guennebaud56686742010-06-24 21:44:24 +0200125 const Matrix saA = A.adjoint() * A;
Gael Guennebaud88e05102012-06-12 13:12:47 +0200126 const Vector b(Vector::Random(size));
127 Vector x(size);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100128
129 // Cholesky module
130 Eigen::LLT<Matrix> LLT; LLT.compute(A);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200131 X = LLT.solve(B);
132 x = LLT.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100133 Eigen::LDLT<Matrix> LDLT; LDLT.compute(A);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200134 X = LDLT.solve(B);
135 x = LDLT.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100136
137 // Eigenvalues module
138 Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp; hessDecomp.compute(complexA);
139 Eigen::ComplexSchur<ComplexMatrix> cSchur(size); cSchur.compute(complexA);
Gael Guennebaud3d8e1792011-01-27 18:02:49 +0100140 Eigen::ComplexEigenSolver<ComplexMatrix> cEigSolver; cEigSolver.compute(complexA);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100141 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);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200147 X = ppLU.solve(B);
148 x = ppLU.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100149 Eigen::FullPivLU<Matrix> fpLU; fpLU.compute(A);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200150 X = fpLU.solve(B);
151 x = fpLU.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100152
153 // QR module
154 Eigen::HouseholderQR<Matrix> hQR; hQR.compute(A);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200155 X = hQR.solve(B);
156 x = hQR.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100157 Eigen::ColPivHouseholderQR<Matrix> cpQR; cpQR.compute(A);
Gael Guennebaud57b58042012-06-26 17:45:01 +0200158 X = cpQR.solve(B);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200159 x = cpQR.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100160 Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A);
Gael Guennebaud882912b2012-06-20 08:58:26 +0200161 // FIXME X = fpQR.solve(B);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200162 x = fpQR.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100163
164 // SVD module
Benoit Jacob8eb0fc12010-10-12 10:19:59 -0400165 Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A, ComputeFullU | ComputeFullV);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100166}
167
Christoph Hertzberg4ba8aa12014-09-25 16:05:17 +0200168void test_zerosized() {
169 // default constructors:
170 Eigen::MatrixXd A;
171 Eigen::VectorXd v;
172 // explicit zero-sized:
173 Eigen::ArrayXXd A0(0,0);
174 Eigen::ArrayXd v0(0);
175
176 // assigning empty objects to each other:
177 A=A0;
178 v=v0;
179}
180
181template<typename MatrixType> void test_reference(const MatrixType& m) {
182 typedef typename MatrixType::Scalar Scalar;
183 enum { Flag = MatrixType::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor};
184 enum { TransposeFlag = !MatrixType::IsRowMajor ? Eigen::RowMajor : Eigen::ColMajor};
185 typename MatrixType::Index rows = m.rows(), cols=m.cols();
Christoph Hertzberg12d59462014-09-30 14:57:54 +0200186 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, Flag > MatrixX;
187 typedef Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic, TransposeFlag> MatrixXT;
Christoph Hertzberg4ba8aa12014-09-25 16:05:17 +0200188 // Dynamic reference:
Christoph Hertzberg12d59462014-09-30 14:57:54 +0200189 typedef Eigen::Ref<const MatrixX > Ref;
190 typedef Eigen::Ref<const MatrixXT > RefT;
Christoph Hertzberg4ba8aa12014-09-25 16:05:17 +0200191
192 Ref r1(m);
193 Ref r2(m.block(rows/3, cols/4, rows/2, cols/2));
194 RefT r3(m.transpose());
195 RefT r4(m.topLeftCorner(rows/2, cols/2).transpose());
196
197 VERIFY_RAISES_ASSERT(RefT r5(m));
198 VERIFY_RAISES_ASSERT(Ref r6(m.transpose()));
199 VERIFY_RAISES_ASSERT(Ref r7(Scalar(2) * m));
Christoph Hertzberg12d59462014-09-30 14:57:54 +0200200
201 // Copy constructors shall also never malloc
202 Ref r8 = r1;
203 RefT r9 = r3;
204
205 // Initializing from a compatible Ref shall also never malloc
206 Eigen::Ref<const MatrixX, Unaligned, Stride<Dynamic, Dynamic> > r10=r8, r11=m;
207
208 // Initializing from an incompatible Ref will malloc:
209 typedef Eigen::Ref<const MatrixX, Aligned> RefAligned;
210 VERIFY_RAISES_ASSERT(RefAligned r12=r10);
211 VERIFY_RAISES_ASSERT(Ref r13=r10); // r10 has more dynamic strides
212
Christoph Hertzberg4ba8aa12014-09-25 16:05:17 +0200213}
214
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +0000215void test_nomalloc()
216{
Christoph Hertzberg12d59462014-09-30 14:57:54 +0200217 // create some dynamic objects
218 Eigen::MatrixXd M1 = MatrixXd::Random(3,3);
219 Ref<const MatrixXd> R1 = 2.0*M1; // Ref requires temporary
220
221 // from here on prohibit malloc:
222 Eigen::internal::set_is_malloc_allowed(false);
223
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +0000224 // check that our operator new is indeed called:
Benoit Jacobaaaade42010-05-30 16:00:58 -0400225 VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3)));
Benoit Jacobb81351c2010-03-08 21:30:06 -0500226 CALL_SUBTEST_1(nomalloc(Matrix<float, 1, 1>()) );
227 CALL_SUBTEST_2(nomalloc(Matrix4d()) );
228 CALL_SUBTEST_3(nomalloc(Matrix<float,32,32>()) );
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100229
230 // Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
Gael Guennebaud8a380472010-07-04 15:35:21 +0200231 CALL_SUBTEST_4(ctms_decompositions<float>());
Christoph Hertzberg12d59462014-09-30 14:57:54 +0200232
Christoph Hertzberg4ba8aa12014-09-25 16:05:17 +0200233 CALL_SUBTEST_5(test_zerosized());
Christoph Hertzberg12d59462014-09-30 14:57:54 +0200234
Christoph Hertzberg4ba8aa12014-09-25 16:05:17 +0200235 CALL_SUBTEST_6(test_reference(Matrix<float,32,32>()));
Christoph Hertzberg12d59462014-09-30 14:57:54 +0200236 CALL_SUBTEST_7(test_reference(R1));
237 CALL_SUBTEST_8(Ref<MatrixXd> R2 = M1.topRows<2>(); test_reference(R2));
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +0000238}