blob: d4ffcefcb1d446dc69496d22d0893a162ddb9907 [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
Benoit Jacob1d52bd42009-01-08 15:20:21 +000015// discard stack allocation as that too bypasses malloc
Benoit Jacobd6712052008-12-31 00:25:31 +000016#define EIGEN_STACK_ALLOCATION_LIMIT 0
Benoit Jacob1d52bd42009-01-08 15:20:21 +000017// any heap allocation will raise an assert
18#define EIGEN_NO_MALLOC
Benoit Jacobd6712052008-12-31 00:25:31 +000019
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000020#include "main.h"
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +010021#include <Eigen/Cholesky>
22#include <Eigen/Eigenvalues>
23#include <Eigen/LU>
24#include <Eigen/QR>
25#include <Eigen/SVD>
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000026
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000027template<typename MatrixType> void nomalloc(const MatrixType& m)
28{
29 /* this test check no dynamic memory allocation are issued with fixed-size matrices
30 */
Hauke Heibelf1679c72010-06-20 17:37:56 +020031 typedef typename MatrixType::Index Index;
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000032 typedef typename MatrixType::Scalar Scalar;
33 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
34
Hauke Heibelf1679c72010-06-20 17:37:56 +020035 Index rows = m.rows();
36 Index cols = m.cols();
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000037
Gael Guennebaudc10f0692008-07-21 00:34:46 +000038 MatrixType m1 = MatrixType::Random(rows, cols),
39 m2 = MatrixType::Random(rows, cols),
Benoit Jacob0609dbe2011-10-31 00:51:36 -040040 m3(rows, cols);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000041
Benoit Jacob47160402010-10-25 10:15:22 -040042 Scalar s1 = internal::random<Scalar>();
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000043
Benoit Jacob47160402010-10-25 10:15:22 -040044 Index r = internal::random<Index>(0, rows-1),
45 c = internal::random<Index>(0, cols-1);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000046
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)));
Gael Guennebaud71a171c2010-01-04 19:13:08 +010049 VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix());
Gael Guennebaud3d8e1792011-01-27 18:02:49 +010050 VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
Gael Guennebaudc60818f2011-02-01 11:38:46 +010051
Gael Guennebaudf46ace62011-01-31 21:30:27 +010052 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();
Gael Guennebaudc60818f2011-02-01 11:38:46 +010061 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);
Gael Guennebaud59af20b2011-02-01 16:46:35 +010084
85 m2.template selfadjointView<Lower>().rankUpdate(m1.col(0),-1);
86 m2.template selfadjointView<Lower>().rankUpdate(m1.row(0),-1);
Gael Guennebaudc60818f2011-02-01 11:38:46 +010087
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);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000094}
95
Gael Guennebaud8a380472010-07-04 15:35:21 +020096template<typename Scalar>
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +010097void ctms_decompositions()
98{
99 const int maxSize = 16;
100 const int size = 12;
101
Gael Guennebaud8a380472010-07-04 15:35:21 +0200102 typedef Eigen::Matrix<Scalar,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100103 Eigen::Dynamic, Eigen::Dynamic,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500104 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100105 maxSize, maxSize> Matrix;
106
Gael Guennebaud8a380472010-07-04 15:35:21 +0200107 typedef Eigen::Matrix<Scalar,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100108 Eigen::Dynamic, 1,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500109 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100110 maxSize, 1> Vector;
111
Gael Guennebaud8a380472010-07-04 15:35:21 +0200112 typedef Eigen::Matrix<std::complex<Scalar>,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100113 Eigen::Dynamic, Eigen::Dynamic,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500114 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100115 maxSize, maxSize> ComplexMatrix;
116
Gael Guennebaud88e05102012-06-12 13:12:47 +0200117 const Matrix A(Matrix::Random(size, size)), B(Matrix::Random(size, size));
118 Matrix X(size,size);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100119 const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
Gael Guennebaud56686742010-06-24 21:44:24 +0200120 const Matrix saA = A.adjoint() * A;
Gael Guennebaud88e05102012-06-12 13:12:47 +0200121 const Vector b(Vector::Random(size));
122 Vector x(size);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100123
124 // Cholesky module
125 Eigen::LLT<Matrix> LLT; LLT.compute(A);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200126 X = LLT.solve(B);
127 x = LLT.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100128 Eigen::LDLT<Matrix> LDLT; LDLT.compute(A);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200129 X = LDLT.solve(B);
130 x = LDLT.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100131
132 // Eigenvalues module
133 Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp; hessDecomp.compute(complexA);
134 Eigen::ComplexSchur<ComplexMatrix> cSchur(size); cSchur.compute(complexA);
Gael Guennebaud3d8e1792011-01-27 18:02:49 +0100135 Eigen::ComplexEigenSolver<ComplexMatrix> cEigSolver; cEigSolver.compute(complexA);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100136 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);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200142 X = ppLU.solve(B);
143 x = ppLU.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100144 Eigen::FullPivLU<Matrix> fpLU; fpLU.compute(A);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200145 X = fpLU.solve(B);
146 x = fpLU.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100147
148 // QR module
149 Eigen::HouseholderQR<Matrix> hQR; hQR.compute(A);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200150 X = hQR.solve(B);
151 x = hQR.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100152 Eigen::ColPivHouseholderQR<Matrix> cpQR; cpQR.compute(A);
Gael Guennebaud57b58042012-06-26 17:45:01 +0200153 X = cpQR.solve(B);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200154 x = cpQR.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100155 Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A);
Gael Guennebaud882912b2012-06-20 08:58:26 +0200156 // FIXME X = fpQR.solve(B);
Gael Guennebaud88e05102012-06-12 13:12:47 +0200157 x = fpQR.solve(b);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100158
159 // SVD module
Benoit Jacob8eb0fc12010-10-12 10:19:59 -0400160 Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A, ComputeFullU | ComputeFullV);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100161}
162
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +0000163void test_nomalloc()
164{
165 // check that our operator new is indeed called:
Benoit Jacobaaaade42010-05-30 16:00:58 -0400166 VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3)));
Benoit Jacobb81351c2010-03-08 21:30:06 -0500167 CALL_SUBTEST_1(nomalloc(Matrix<float, 1, 1>()) );
168 CALL_SUBTEST_2(nomalloc(Matrix4d()) );
169 CALL_SUBTEST_3(nomalloc(Matrix<float,32,32>()) );
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100170
171 // Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
Gael Guennebaud8a380472010-07-04 15:35:21 +0200172 CALL_SUBTEST_4(ctms_decompositions<float>());
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100173
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +0000174}