2008-06-02 22:54:52 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-23 02:25:33 +08:00
|
|
|
// for linear algebra.
|
2008-06-02 22:54:52 +08:00
|
|
|
//
|
2010-06-25 05:21:58 +08:00
|
|
|
// Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
|
2008-11-24 21:40:43 +08:00
|
|
|
// Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
|
2008-06-02 22:54:52 +08:00
|
|
|
//
|
2012-07-14 02:42:47 +08:00
|
|
|
// This Source Code Form is subject to the terms of the Mozilla
|
|
|
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
|
|
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
2008-06-02 22:54:52 +08:00
|
|
|
|
2008-06-08 23:03:23 +08:00
|
|
|
// this hack is needed to make this file compiles with -pedantic (gcc)
|
2009-01-22 08:09:34 +08:00
|
|
|
#ifdef __GNUC__
|
2008-06-04 18:15:48 +08:00
|
|
|
#define throw(X)
|
2009-01-22 08:09:34 +08:00
|
|
|
#endif
|
2009-01-08 23:20:21 +08:00
|
|
|
// discard stack allocation as that too bypasses malloc
|
2008-12-31 08:25:31 +08:00
|
|
|
#define EIGEN_STACK_ALLOCATION_LIMIT 0
|
2009-01-08 23:20:21 +08:00
|
|
|
// any heap allocation will raise an assert
|
|
|
|
#define EIGEN_NO_MALLOC
|
2008-12-31 08:25:31 +08:00
|
|
|
|
2008-06-02 22:54:52 +08:00
|
|
|
#include "main.h"
|
2010-03-09 02:31:27 +08:00
|
|
|
#include <Eigen/Cholesky>
|
|
|
|
#include <Eigen/Eigenvalues>
|
|
|
|
#include <Eigen/LU>
|
|
|
|
#include <Eigen/QR>
|
|
|
|
#include <Eigen/SVD>
|
2008-06-02 22:54:52 +08:00
|
|
|
|
|
|
|
template<typename MatrixType> void nomalloc(const MatrixType& m)
|
|
|
|
{
|
|
|
|
/* this test check no dynamic memory allocation are issued with fixed-size matrices
|
|
|
|
*/
|
2010-06-20 23:37:56 +08:00
|
|
|
typedef typename MatrixType::Index Index;
|
2008-06-02 22:54:52 +08:00
|
|
|
typedef typename MatrixType::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
|
|
|
|
|
2010-06-20 23:37:56 +08:00
|
|
|
Index rows = m.rows();
|
|
|
|
Index cols = m.cols();
|
2008-06-02 22:54:52 +08:00
|
|
|
|
2008-07-21 08:34:46 +08:00
|
|
|
MatrixType m1 = MatrixType::Random(rows, cols),
|
|
|
|
m2 = MatrixType::Random(rows, cols),
|
2011-10-31 12:51:36 +08:00
|
|
|
m3(rows, cols);
|
2008-06-02 22:54:52 +08:00
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
Scalar s1 = internal::random<Scalar>();
|
2008-06-02 22:54:52 +08:00
|
|
|
|
2010-10-25 22:15:22 +08:00
|
|
|
Index r = internal::random<Index>(0, rows-1),
|
|
|
|
c = internal::random<Index>(0, cols-1);
|
2008-06-02 22:54:52 +08:00
|
|
|
|
|
|
|
VERIFY_IS_APPROX((m1+m2)*s1, s1*m1+s1*m2);
|
|
|
|
VERIFY_IS_APPROX((m1+m2)(r,c), (m1(r,c))+(m2(r,c)));
|
2010-01-05 02:13:08 +08:00
|
|
|
VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix());
|
2011-01-28 01:02:49 +08:00
|
|
|
VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
|
2011-02-01 18:38:46 +08:00
|
|
|
|
2011-02-01 04:30:27 +08:00
|
|
|
m2.col(0).noalias() = m1 * m1.col(0);
|
|
|
|
m2.col(0).noalias() -= m1.adjoint() * m1.col(0);
|
|
|
|
m2.col(0).noalias() -= m1 * m1.row(0).adjoint();
|
|
|
|
m2.col(0).noalias() -= m1.adjoint() * m1.row(0).adjoint();
|
|
|
|
|
|
|
|
m2.row(0).noalias() = m1.row(0) * m1;
|
|
|
|
m2.row(0).noalias() -= m1.row(0) * m1.adjoint();
|
|
|
|
m2.row(0).noalias() -= m1.col(0).adjoint() * m1;
|
|
|
|
m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint();
|
2011-02-01 18:38:46 +08:00
|
|
|
VERIFY_IS_APPROX(m2,m2);
|
|
|
|
|
|
|
|
m2.col(0).noalias() = m1.template triangularView<Upper>() * m1.col(0);
|
|
|
|
m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.col(0);
|
|
|
|
m2.col(0).noalias() -= m1.template triangularView<Upper>() * m1.row(0).adjoint();
|
|
|
|
m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.row(0).adjoint();
|
|
|
|
|
|
|
|
m2.row(0).noalias() = m1.row(0) * m1.template triangularView<Upper>();
|
|
|
|
m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template triangularView<Upper>();
|
|
|
|
m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template triangularView<Upper>();
|
|
|
|
m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template triangularView<Upper>();
|
|
|
|
VERIFY_IS_APPROX(m2,m2);
|
|
|
|
|
|
|
|
m2.col(0).noalias() = m1.template selfadjointView<Upper>() * m1.col(0);
|
|
|
|
m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.col(0);
|
|
|
|
m2.col(0).noalias() -= m1.template selfadjointView<Upper>() * m1.row(0).adjoint();
|
|
|
|
m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.row(0).adjoint();
|
|
|
|
|
|
|
|
m2.row(0).noalias() = m1.row(0) * m1.template selfadjointView<Upper>();
|
|
|
|
m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template selfadjointView<Upper>();
|
|
|
|
m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template selfadjointView<Upper>();
|
|
|
|
m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template selfadjointView<Upper>();
|
|
|
|
VERIFY_IS_APPROX(m2,m2);
|
2011-02-01 23:46:35 +08:00
|
|
|
|
|
|
|
m2.template selfadjointView<Lower>().rankUpdate(m1.col(0),-1);
|
|
|
|
m2.template selfadjointView<Lower>().rankUpdate(m1.row(0),-1);
|
2011-02-01 18:38:46 +08:00
|
|
|
|
|
|
|
// The following fancy matrix-matrix products are not safe yet regarding static allocation
|
|
|
|
// m1 += m1.template triangularView<Upper>() * m2.col(;
|
|
|
|
// m1.template selfadjointView<Lower>().rankUpdate(m2);
|
|
|
|
// m1 += m1.template triangularView<Upper>() * m2;
|
|
|
|
// m1 += m1.template selfadjointView<Lower>() * m2;
|
|
|
|
// VERIFY_IS_APPROX(m1,m1);
|
2008-06-02 22:54:52 +08:00
|
|
|
}
|
|
|
|
|
2010-07-04 21:35:21 +08:00
|
|
|
template<typename Scalar>
|
2010-03-09 02:31:27 +08:00
|
|
|
void ctms_decompositions()
|
|
|
|
{
|
|
|
|
const int maxSize = 16;
|
|
|
|
const int size = 12;
|
|
|
|
|
2010-07-04 21:35:21 +08:00
|
|
|
typedef Eigen::Matrix<Scalar,
|
2010-03-09 02:31:27 +08:00
|
|
|
Eigen::Dynamic, Eigen::Dynamic,
|
2010-03-09 10:30:06 +08:00
|
|
|
0,
|
2010-03-09 02:31:27 +08:00
|
|
|
maxSize, maxSize> Matrix;
|
|
|
|
|
2010-07-04 21:35:21 +08:00
|
|
|
typedef Eigen::Matrix<Scalar,
|
2010-03-09 02:31:27 +08:00
|
|
|
Eigen::Dynamic, 1,
|
2010-03-09 10:30:06 +08:00
|
|
|
0,
|
2010-03-09 02:31:27 +08:00
|
|
|
maxSize, 1> Vector;
|
|
|
|
|
2010-07-04 21:35:21 +08:00
|
|
|
typedef Eigen::Matrix<std::complex<Scalar>,
|
2010-03-09 02:31:27 +08:00
|
|
|
Eigen::Dynamic, Eigen::Dynamic,
|
2010-03-09 10:30:06 +08:00
|
|
|
0,
|
2010-03-09 02:31:27 +08:00
|
|
|
maxSize, maxSize> ComplexMatrix;
|
|
|
|
|
2012-06-12 19:12:47 +08:00
|
|
|
const Matrix A(Matrix::Random(size, size)), B(Matrix::Random(size, size));
|
|
|
|
Matrix X(size,size);
|
2010-03-09 02:31:27 +08:00
|
|
|
const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
|
2010-06-25 03:44:24 +08:00
|
|
|
const Matrix saA = A.adjoint() * A;
|
2012-06-12 19:12:47 +08:00
|
|
|
const Vector b(Vector::Random(size));
|
|
|
|
Vector x(size);
|
2010-03-09 02:31:27 +08:00
|
|
|
|
|
|
|
// Cholesky module
|
|
|
|
Eigen::LLT<Matrix> LLT; LLT.compute(A);
|
2012-06-12 19:12:47 +08:00
|
|
|
X = LLT.solve(B);
|
|
|
|
x = LLT.solve(b);
|
2010-03-09 02:31:27 +08:00
|
|
|
Eigen::LDLT<Matrix> LDLT; LDLT.compute(A);
|
2012-06-12 19:12:47 +08:00
|
|
|
X = LDLT.solve(B);
|
|
|
|
x = LDLT.solve(b);
|
2010-03-09 02:31:27 +08:00
|
|
|
|
|
|
|
// Eigenvalues module
|
|
|
|
Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp; hessDecomp.compute(complexA);
|
|
|
|
Eigen::ComplexSchur<ComplexMatrix> cSchur(size); cSchur.compute(complexA);
|
2011-01-28 01:02:49 +08:00
|
|
|
Eigen::ComplexEigenSolver<ComplexMatrix> cEigSolver; cEigSolver.compute(complexA);
|
2010-03-09 02:31:27 +08:00
|
|
|
Eigen::EigenSolver<Matrix> eigSolver; eigSolver.compute(A);
|
|
|
|
Eigen::SelfAdjointEigenSolver<Matrix> saEigSolver(size); saEigSolver.compute(saA);
|
|
|
|
Eigen::Tridiagonalization<Matrix> tridiag; tridiag.compute(saA);
|
|
|
|
|
|
|
|
// LU module
|
|
|
|
Eigen::PartialPivLU<Matrix> ppLU; ppLU.compute(A);
|
2012-06-12 19:12:47 +08:00
|
|
|
X = ppLU.solve(B);
|
|
|
|
x = ppLU.solve(b);
|
2010-03-09 02:31:27 +08:00
|
|
|
Eigen::FullPivLU<Matrix> fpLU; fpLU.compute(A);
|
2012-06-12 19:12:47 +08:00
|
|
|
X = fpLU.solve(B);
|
|
|
|
x = fpLU.solve(b);
|
2010-03-09 02:31:27 +08:00
|
|
|
|
|
|
|
// QR module
|
|
|
|
Eigen::HouseholderQR<Matrix> hQR; hQR.compute(A);
|
2012-06-12 19:12:47 +08:00
|
|
|
X = hQR.solve(B);
|
|
|
|
x = hQR.solve(b);
|
2010-03-09 02:31:27 +08:00
|
|
|
Eigen::ColPivHouseholderQR<Matrix> cpQR; cpQR.compute(A);
|
2012-06-26 23:45:01 +08:00
|
|
|
X = cpQR.solve(B);
|
2012-06-12 19:12:47 +08:00
|
|
|
x = cpQR.solve(b);
|
2010-03-09 02:31:27 +08:00
|
|
|
Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A);
|
2012-06-20 14:58:26 +08:00
|
|
|
// FIXME X = fpQR.solve(B);
|
2012-06-12 19:12:47 +08:00
|
|
|
x = fpQR.solve(b);
|
2010-03-09 02:31:27 +08:00
|
|
|
|
|
|
|
// SVD module
|
2010-10-12 22:19:59 +08:00
|
|
|
Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A, ComputeFullU | ComputeFullV);
|
2010-03-09 02:31:27 +08:00
|
|
|
}
|
|
|
|
|
2008-06-02 22:54:52 +08:00
|
|
|
void test_nomalloc()
|
|
|
|
{
|
|
|
|
// check that our operator new is indeed called:
|
2010-05-31 04:00:58 +08:00
|
|
|
VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3)));
|
2010-03-09 10:30:06 +08:00
|
|
|
CALL_SUBTEST_1(nomalloc(Matrix<float, 1, 1>()) );
|
|
|
|
CALL_SUBTEST_2(nomalloc(Matrix4d()) );
|
|
|
|
CALL_SUBTEST_3(nomalloc(Matrix<float,32,32>()) );
|
2010-03-09 02:31:27 +08:00
|
|
|
|
|
|
|
// Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
|
2010-07-04 21:35:21 +08:00
|
|
|
CALL_SUBTEST_4(ctms_decompositions<float>());
|
2010-03-09 02:31:27 +08:00
|
|
|
|
2008-06-02 22:54:52 +08:00
|
|
|
}
|