blob: 7ef71bfcd3440e2b5dec220b7075b416bec62c46 [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//
7// Eigen is free software; you can redistribute it and/or
8// modify it under the terms of the GNU Lesser General Public
9// License as published by the Free Software Foundation; either
10// version 3 of the License, or (at your option) any later version.
11//
12// Alternatively, you can redistribute it and/or
13// modify it under the terms of the GNU General Public License as
14// published by the Free Software Foundation; either version 2 of
15// the License, or (at your option) any later version.
16//
17// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20// GNU General Public License for more details.
21//
22// You should have received a copy of the GNU Lesser General Public
23// License and a copy of the GNU General Public License along with
24// Eigen. If not, see <http://www.gnu.org/licenses/>.
25
Gael Guennebaude3fac692008-06-08 15:03:23 +000026// this hack is needed to make this file compiles with -pedantic (gcc)
Benoit Jacob9e3c7312009-01-22 00:09:34 +000027#ifdef __GNUC__
Gael Guennebaud60726f92008-06-04 10:15:48 +000028#define throw(X)
Benoit Jacob9e3c7312009-01-22 00:09:34 +000029#endif
Benoit Jacob1d52bd42009-01-08 15:20:21 +000030// discard stack allocation as that too bypasses malloc
Benoit Jacobd6712052008-12-31 00:25:31 +000031#define EIGEN_STACK_ALLOCATION_LIMIT 0
Benoit Jacob1d52bd42009-01-08 15:20:21 +000032// any heap allocation will raise an assert
33#define EIGEN_NO_MALLOC
Benoit Jacobd6712052008-12-31 00:25:31 +000034
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000035#include "main.h"
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +010036#include <Eigen/Cholesky>
37#include <Eigen/Eigenvalues>
38#include <Eigen/LU>
39#include <Eigen/QR>
40#include <Eigen/SVD>
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000041
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000042template<typename MatrixType> void nomalloc(const MatrixType& m)
43{
44 /* this test check no dynamic memory allocation are issued with fixed-size matrices
45 */
Hauke Heibelf1679c72010-06-20 17:37:56 +020046 typedef typename MatrixType::Index Index;
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000047 typedef typename MatrixType::Scalar Scalar;
48 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
49
Hauke Heibelf1679c72010-06-20 17:37:56 +020050 Index rows = m.rows();
51 Index cols = m.cols();
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000052
Gael Guennebaudc10f0692008-07-21 00:34:46 +000053 MatrixType m1 = MatrixType::Random(rows, cols),
54 m2 = MatrixType::Random(rows, cols),
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000055 m3(rows, cols),
Gael Guennebaudc10f0692008-07-21 00:34:46 +000056 mzero = MatrixType::Zero(rows, cols),
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000057 identity = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
Gael Guennebaudc10f0692008-07-21 00:34:46 +000058 ::Identity(rows, rows),
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000059 square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>
Gael Guennebaudc10f0692008-07-21 00:34:46 +000060 ::Random(rows, rows);
61 VectorType v1 = VectorType::Random(rows),
62 v2 = VectorType::Random(rows),
63 vzero = VectorType::Zero(rows);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000064
Benoit Jacob47160402010-10-25 10:15:22 -040065 Scalar s1 = internal::random<Scalar>();
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000066
Benoit Jacob47160402010-10-25 10:15:22 -040067 Index r = internal::random<Index>(0, rows-1),
68 c = internal::random<Index>(0, cols-1);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +000069
70 VERIFY_IS_APPROX((m1+m2)*s1, s1*m1+s1*m2);
71 VERIFY_IS_APPROX((m1+m2)(r,c), (m1(r,c))+(m2(r,c)));
Gael Guennebaud71a171c2010-01-04 19:13:08 +010072 VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix());
Gael Guennebaud3d8e1792011-01-27 18:02:49 +010073 VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
Gael Guennebaudc60818f2011-02-01 11:38:46 +010074
Gael Guennebaudf46ace62011-01-31 21:30:27 +010075 m2.col(0).noalias() = m1 * m1.col(0);
76 m2.col(0).noalias() -= m1.adjoint() * m1.col(0);
77 m2.col(0).noalias() -= m1 * m1.row(0).adjoint();
78 m2.col(0).noalias() -= m1.adjoint() * m1.row(0).adjoint();
79
80 m2.row(0).noalias() = m1.row(0) * m1;
81 m2.row(0).noalias() -= m1.row(0) * m1.adjoint();
82 m2.row(0).noalias() -= m1.col(0).adjoint() * m1;
83 m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint();
Gael Guennebaudc60818f2011-02-01 11:38:46 +010084 VERIFY_IS_APPROX(m2,m2);
85
86 m2.col(0).noalias() = m1.template triangularView<Upper>() * m1.col(0);
87 m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.col(0);
88 m2.col(0).noalias() -= m1.template triangularView<Upper>() * m1.row(0).adjoint();
89 m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.row(0).adjoint();
90
91 m2.row(0).noalias() = m1.row(0) * m1.template triangularView<Upper>();
92 m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template triangularView<Upper>();
93 m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template triangularView<Upper>();
94 m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template triangularView<Upper>();
95 VERIFY_IS_APPROX(m2,m2);
96
97 m2.col(0).noalias() = m1.template selfadjointView<Upper>() * m1.col(0);
98 m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.col(0);
99 m2.col(0).noalias() -= m1.template selfadjointView<Upper>() * m1.row(0).adjoint();
100 m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.row(0).adjoint();
101
102 m2.row(0).noalias() = m1.row(0) * m1.template selfadjointView<Upper>();
103 m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template selfadjointView<Upper>();
104 m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template selfadjointView<Upper>();
105 m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template selfadjointView<Upper>();
106 VERIFY_IS_APPROX(m2,m2);
107
108 // The following fancy matrix-matrix products are not safe yet regarding static allocation
109// m1 += m1.template triangularView<Upper>() * m2.col(;
110// m1.template selfadjointView<Lower>().rankUpdate(m2);
111// m1 += m1.template triangularView<Upper>() * m2;
112// m1 += m1.template selfadjointView<Lower>() * m2;
113// VERIFY_IS_APPROX(m1,m1);
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +0000114}
115
Gael Guennebaud8a380472010-07-04 15:35:21 +0200116template<typename Scalar>
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100117void ctms_decompositions()
118{
119 const int maxSize = 16;
120 const int size = 12;
121
Gael Guennebaud8a380472010-07-04 15:35:21 +0200122 typedef Eigen::Matrix<Scalar,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100123 Eigen::Dynamic, Eigen::Dynamic,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500124 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100125 maxSize, maxSize> Matrix;
126
Gael Guennebaud8a380472010-07-04 15:35:21 +0200127 typedef Eigen::Matrix<Scalar,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100128 Eigen::Dynamic, 1,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500129 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100130 maxSize, 1> Vector;
131
Gael Guennebaud8a380472010-07-04 15:35:21 +0200132 typedef Eigen::Matrix<std::complex<Scalar>,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100133 Eigen::Dynamic, Eigen::Dynamic,
Benoit Jacobb81351c2010-03-08 21:30:06 -0500134 0,
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100135 maxSize, maxSize> ComplexMatrix;
136
137 const Matrix A(Matrix::Random(size, size));
138 const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
Gael Guennebaud56686742010-06-24 21:44:24 +0200139 const Matrix saA = A.adjoint() * A;
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100140
141 // Cholesky module
142 Eigen::LLT<Matrix> LLT; LLT.compute(A);
143 Eigen::LDLT<Matrix> LDLT; LDLT.compute(A);
144
145 // Eigenvalues module
146 Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp; hessDecomp.compute(complexA);
147 Eigen::ComplexSchur<ComplexMatrix> cSchur(size); cSchur.compute(complexA);
Gael Guennebaud3d8e1792011-01-27 18:02:49 +0100148 Eigen::ComplexEigenSolver<ComplexMatrix> cEigSolver; cEigSolver.compute(complexA);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100149 Eigen::EigenSolver<Matrix> eigSolver; eigSolver.compute(A);
150 Eigen::SelfAdjointEigenSolver<Matrix> saEigSolver(size); saEigSolver.compute(saA);
151 Eigen::Tridiagonalization<Matrix> tridiag; tridiag.compute(saA);
152
153 // LU module
154 Eigen::PartialPivLU<Matrix> ppLU; ppLU.compute(A);
155 Eigen::FullPivLU<Matrix> fpLU; fpLU.compute(A);
156
157 // QR module
158 Eigen::HouseholderQR<Matrix> hQR; hQR.compute(A);
159 Eigen::ColPivHouseholderQR<Matrix> cpQR; cpQR.compute(A);
160 Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A);
161
162 // SVD module
Benoit Jacob8eb0fc12010-10-12 10:19:59 -0400163 Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A, ComputeFullU | ComputeFullV);
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100164}
165
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +0000166void test_nomalloc()
167{
168 // check that our operator new is indeed called:
Benoit Jacobaaaade42010-05-30 16:00:58 -0400169 VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3)));
Benoit Jacobb81351c2010-03-08 21:30:06 -0500170 CALL_SUBTEST_1(nomalloc(Matrix<float, 1, 1>()) );
171 CALL_SUBTEST_2(nomalloc(Matrix4d()) );
172 CALL_SUBTEST_3(nomalloc(Matrix<float,32,32>()) );
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100173
174 // Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
Gael Guennebaud8a380472010-07-04 15:35:21 +0200175 CALL_SUBTEST_4(ctms_decompositions<float>());
Adolfo Rodriguez Tsouroukdissian5a36f4a2010-03-08 19:31:27 +0100176
Gael Guennebaudf2ebbee2008-06-02 14:54:52 +0000177}