Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 1 | // This file is part of Eigen, a lightweight C++ template library |
Benoit Jacob | 6347b1d | 2009-05-22 20:25:33 +0200 | [diff] [blame] | 2 | // for linear algebra. |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 3 | // |
Gael Guennebaud | 28e64b0 | 2010-06-24 23:21:58 +0200 | [diff] [blame] | 4 | // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr> |
Benoit Jacob | 00f89a8 | 2008-11-24 13:40:43 +0000 | [diff] [blame] | 5 | // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com> |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 6 | // |
Benoit Jacob | 69124cf | 2012-07-13 14:42:47 -0400 | [diff] [blame^] | 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/. |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 10 | |
Gael Guennebaud | e3fac69 | 2008-06-08 15:03:23 +0000 | [diff] [blame] | 11 | // this hack is needed to make this file compiles with -pedantic (gcc) |
Benoit Jacob | 9e3c731 | 2009-01-22 00:09:34 +0000 | [diff] [blame] | 12 | #ifdef __GNUC__ |
Gael Guennebaud | 60726f9 | 2008-06-04 10:15:48 +0000 | [diff] [blame] | 13 | #define throw(X) |
Benoit Jacob | 9e3c731 | 2009-01-22 00:09:34 +0000 | [diff] [blame] | 14 | #endif |
Benoit Jacob | 1d52bd4 | 2009-01-08 15:20:21 +0000 | [diff] [blame] | 15 | // discard stack allocation as that too bypasses malloc |
Benoit Jacob | d671205 | 2008-12-31 00:25:31 +0000 | [diff] [blame] | 16 | #define EIGEN_STACK_ALLOCATION_LIMIT 0 |
Benoit Jacob | 1d52bd4 | 2009-01-08 15:20:21 +0000 | [diff] [blame] | 17 | // any heap allocation will raise an assert |
| 18 | #define EIGEN_NO_MALLOC |
Benoit Jacob | d671205 | 2008-12-31 00:25:31 +0000 | [diff] [blame] | 19 | |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 20 | #include "main.h" |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 21 | #include <Eigen/Cholesky> |
| 22 | #include <Eigen/Eigenvalues> |
| 23 | #include <Eigen/LU> |
| 24 | #include <Eigen/QR> |
| 25 | #include <Eigen/SVD> |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 26 | |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 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 | */ |
Hauke Heibel | f1679c7 | 2010-06-20 17:37:56 +0200 | [diff] [blame] | 31 | typedef typename MatrixType::Index Index; |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 32 | typedef typename MatrixType::Scalar Scalar; |
| 33 | typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType; |
| 34 | |
Hauke Heibel | f1679c7 | 2010-06-20 17:37:56 +0200 | [diff] [blame] | 35 | Index rows = m.rows(); |
| 36 | Index cols = m.cols(); |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 37 | |
Gael Guennebaud | c10f069 | 2008-07-21 00:34:46 +0000 | [diff] [blame] | 38 | MatrixType m1 = MatrixType::Random(rows, cols), |
| 39 | m2 = MatrixType::Random(rows, cols), |
Benoit Jacob | 0609dbe | 2011-10-31 00:51:36 -0400 | [diff] [blame] | 40 | m3(rows, cols); |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 41 | |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 42 | Scalar s1 = internal::random<Scalar>(); |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 43 | |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 44 | Index r = internal::random<Index>(0, rows-1), |
| 45 | c = internal::random<Index>(0, cols-1); |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 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))); |
Gael Guennebaud | 71a171c | 2010-01-04 19:13:08 +0100 | [diff] [blame] | 49 | VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix()); |
Gael Guennebaud | 3d8e179 | 2011-01-27 18:02:49 +0100 | [diff] [blame] | 50 | VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2)); |
Gael Guennebaud | c60818f | 2011-02-01 11:38:46 +0100 | [diff] [blame] | 51 | |
Gael Guennebaud | f46ace6 | 2011-01-31 21:30:27 +0100 | [diff] [blame] | 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(); |
Gael Guennebaud | c60818f | 2011-02-01 11:38:46 +0100 | [diff] [blame] | 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); |
Gael Guennebaud | 59af20b | 2011-02-01 16:46:35 +0100 | [diff] [blame] | 84 | |
| 85 | m2.template selfadjointView<Lower>().rankUpdate(m1.col(0),-1); |
| 86 | m2.template selfadjointView<Lower>().rankUpdate(m1.row(0),-1); |
Gael Guennebaud | c60818f | 2011-02-01 11:38:46 +0100 | [diff] [blame] | 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); |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 94 | } |
| 95 | |
Gael Guennebaud | 8a38047 | 2010-07-04 15:35:21 +0200 | [diff] [blame] | 96 | template<typename Scalar> |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 97 | void ctms_decompositions() |
| 98 | { |
| 99 | const int maxSize = 16; |
| 100 | const int size = 12; |
| 101 | |
Gael Guennebaud | 8a38047 | 2010-07-04 15:35:21 +0200 | [diff] [blame] | 102 | typedef Eigen::Matrix<Scalar, |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 103 | Eigen::Dynamic, Eigen::Dynamic, |
Benoit Jacob | b81351c | 2010-03-08 21:30:06 -0500 | [diff] [blame] | 104 | 0, |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 105 | maxSize, maxSize> Matrix; |
| 106 | |
Gael Guennebaud | 8a38047 | 2010-07-04 15:35:21 +0200 | [diff] [blame] | 107 | typedef Eigen::Matrix<Scalar, |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 108 | Eigen::Dynamic, 1, |
Benoit Jacob | b81351c | 2010-03-08 21:30:06 -0500 | [diff] [blame] | 109 | 0, |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 110 | maxSize, 1> Vector; |
| 111 | |
Gael Guennebaud | 8a38047 | 2010-07-04 15:35:21 +0200 | [diff] [blame] | 112 | typedef Eigen::Matrix<std::complex<Scalar>, |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 113 | Eigen::Dynamic, Eigen::Dynamic, |
Benoit Jacob | b81351c | 2010-03-08 21:30:06 -0500 | [diff] [blame] | 114 | 0, |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 115 | maxSize, maxSize> ComplexMatrix; |
| 116 | |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 117 | const Matrix A(Matrix::Random(size, size)), B(Matrix::Random(size, size)); |
| 118 | Matrix X(size,size); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 119 | const ComplexMatrix complexA(ComplexMatrix::Random(size, size)); |
Gael Guennebaud | 5668674 | 2010-06-24 21:44:24 +0200 | [diff] [blame] | 120 | const Matrix saA = A.adjoint() * A; |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 121 | const Vector b(Vector::Random(size)); |
| 122 | Vector x(size); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 123 | |
| 124 | // Cholesky module |
| 125 | Eigen::LLT<Matrix> LLT; LLT.compute(A); |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 126 | X = LLT.solve(B); |
| 127 | x = LLT.solve(b); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 128 | Eigen::LDLT<Matrix> LDLT; LDLT.compute(A); |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 129 | X = LDLT.solve(B); |
| 130 | x = LDLT.solve(b); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 131 | |
| 132 | // Eigenvalues module |
| 133 | Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp; hessDecomp.compute(complexA); |
| 134 | Eigen::ComplexSchur<ComplexMatrix> cSchur(size); cSchur.compute(complexA); |
Gael Guennebaud | 3d8e179 | 2011-01-27 18:02:49 +0100 | [diff] [blame] | 135 | Eigen::ComplexEigenSolver<ComplexMatrix> cEigSolver; cEigSolver.compute(complexA); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 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); |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 142 | X = ppLU.solve(B); |
| 143 | x = ppLU.solve(b); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 144 | Eigen::FullPivLU<Matrix> fpLU; fpLU.compute(A); |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 145 | X = fpLU.solve(B); |
| 146 | x = fpLU.solve(b); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 147 | |
| 148 | // QR module |
| 149 | Eigen::HouseholderQR<Matrix> hQR; hQR.compute(A); |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 150 | X = hQR.solve(B); |
| 151 | x = hQR.solve(b); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 152 | Eigen::ColPivHouseholderQR<Matrix> cpQR; cpQR.compute(A); |
Gael Guennebaud | 57b5804 | 2012-06-26 17:45:01 +0200 | [diff] [blame] | 153 | X = cpQR.solve(B); |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 154 | x = cpQR.solve(b); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 155 | Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A); |
Gael Guennebaud | 882912b | 2012-06-20 08:58:26 +0200 | [diff] [blame] | 156 | // FIXME X = fpQR.solve(B); |
Gael Guennebaud | 88e0510 | 2012-06-12 13:12:47 +0200 | [diff] [blame] | 157 | x = fpQR.solve(b); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 158 | |
| 159 | // SVD module |
Benoit Jacob | 8eb0fc1 | 2010-10-12 10:19:59 -0400 | [diff] [blame] | 160 | Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A, ComputeFullU | ComputeFullV); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 161 | } |
| 162 | |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 163 | void test_nomalloc() |
| 164 | { |
| 165 | // check that our operator new is indeed called: |
Benoit Jacob | aaaade4 | 2010-05-30 16:00:58 -0400 | [diff] [blame] | 166 | VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3))); |
Benoit Jacob | b81351c | 2010-03-08 21:30:06 -0500 | [diff] [blame] | 167 | 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 Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 170 | |
| 171 | // Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms) |
Gael Guennebaud | 8a38047 | 2010-07-04 15:35:21 +0200 | [diff] [blame] | 172 | CALL_SUBTEST_4(ctms_decompositions<float>()); |
Adolfo Rodriguez Tsouroukdissian | 5a36f4a | 2010-03-08 19:31:27 +0100 | [diff] [blame] | 173 | |
Gael Guennebaud | f2ebbee | 2008-06-02 14:54:52 +0000 | [diff] [blame] | 174 | } |