You can not select more than 25 topics Topics must start with a letter or number, can include dashes ('-') and can be up to 35 characters long.

174 lines
6.6 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
  6. //
  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/.
  10. // this hack is needed to make this file compiles with -pedantic (gcc)
  11. #ifdef __GNUC__
  12. #define throw(X)
  13. #endif
  14. // discard stack allocation as that too bypasses malloc
  15. #define EIGEN_STACK_ALLOCATION_LIMIT 0
  16. // any heap allocation will raise an assert
  17. #define EIGEN_NO_MALLOC
  18. #include "main.h"
  19. #include <Eigen/Cholesky>
  20. #include <Eigen/Eigenvalues>
  21. #include <Eigen/LU>
  22. #include <Eigen/QR>
  23. #include <Eigen/SVD>
  24. template<typename MatrixType> void nomalloc(const MatrixType& m)
  25. {
  26. /* this test check no dynamic memory allocation are issued with fixed-size matrices
  27. */
  28. typedef typename MatrixType::Index Index;
  29. typedef typename MatrixType::Scalar Scalar;
  30. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  31. Index rows = m.rows();
  32. Index cols = m.cols();
  33. MatrixType m1 = MatrixType::Random(rows, cols),
  34. m2 = MatrixType::Random(rows, cols),
  35. m3(rows, cols);
  36. Scalar s1 = internal::random<Scalar>();
  37. Index r = internal::random<Index>(0, rows-1),
  38. c = internal::random<Index>(0, cols-1);
  39. VERIFY_IS_APPROX((m1+m2)*s1, s1*m1+s1*m2);
  40. VERIFY_IS_APPROX((m1+m2)(r,c), (m1(r,c))+(m2(r,c)));
  41. VERIFY_IS_APPROX(m1.cwiseProduct(m1.block(0,0,rows,cols)), (m1.array()*m1.array()).matrix());
  42. VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
  43. m2.col(0).noalias() = m1 * m1.col(0);
  44. m2.col(0).noalias() -= m1.adjoint() * m1.col(0);
  45. m2.col(0).noalias() -= m1 * m1.row(0).adjoint();
  46. m2.col(0).noalias() -= m1.adjoint() * m1.row(0).adjoint();
  47. m2.row(0).noalias() = m1.row(0) * m1;
  48. m2.row(0).noalias() -= m1.row(0) * m1.adjoint();
  49. m2.row(0).noalias() -= m1.col(0).adjoint() * m1;
  50. m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint();
  51. VERIFY_IS_APPROX(m2,m2);
  52. m2.col(0).noalias() = m1.template triangularView<Upper>() * m1.col(0);
  53. m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.col(0);
  54. m2.col(0).noalias() -= m1.template triangularView<Upper>() * m1.row(0).adjoint();
  55. m2.col(0).noalias() -= m1.adjoint().template triangularView<Upper>() * m1.row(0).adjoint();
  56. m2.row(0).noalias() = m1.row(0) * m1.template triangularView<Upper>();
  57. m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template triangularView<Upper>();
  58. m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template triangularView<Upper>();
  59. m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template triangularView<Upper>();
  60. VERIFY_IS_APPROX(m2,m2);
  61. m2.col(0).noalias() = m1.template selfadjointView<Upper>() * m1.col(0);
  62. m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.col(0);
  63. m2.col(0).noalias() -= m1.template selfadjointView<Upper>() * m1.row(0).adjoint();
  64. m2.col(0).noalias() -= m1.adjoint().template selfadjointView<Upper>() * m1.row(0).adjoint();
  65. m2.row(0).noalias() = m1.row(0) * m1.template selfadjointView<Upper>();
  66. m2.row(0).noalias() -= m1.row(0) * m1.adjoint().template selfadjointView<Upper>();
  67. m2.row(0).noalias() -= m1.col(0).adjoint() * m1.template selfadjointView<Upper>();
  68. m2.row(0).noalias() -= m1.col(0).adjoint() * m1.adjoint().template selfadjointView<Upper>();
  69. VERIFY_IS_APPROX(m2,m2);
  70. m2.template selfadjointView<Lower>().rankUpdate(m1.col(0),-1);
  71. m2.template selfadjointView<Lower>().rankUpdate(m1.row(0),-1);
  72. // The following fancy matrix-matrix products are not safe yet regarding static allocation
  73. // m1 += m1.template triangularView<Upper>() * m2.col(;
  74. // m1.template selfadjointView<Lower>().rankUpdate(m2);
  75. // m1 += m1.template triangularView<Upper>() * m2;
  76. // m1 += m1.template selfadjointView<Lower>() * m2;
  77. // VERIFY_IS_APPROX(m1,m1);
  78. }
  79. template<typename Scalar>
  80. void ctms_decompositions()
  81. {
  82. const int maxSize = 16;
  83. const int size = 12;
  84. typedef Eigen::Matrix<Scalar,
  85. Eigen::Dynamic, Eigen::Dynamic,
  86. 0,
  87. maxSize, maxSize> Matrix;
  88. typedef Eigen::Matrix<Scalar,
  89. Eigen::Dynamic, 1,
  90. 0,
  91. maxSize, 1> Vector;
  92. typedef Eigen::Matrix<std::complex<Scalar>,
  93. Eigen::Dynamic, Eigen::Dynamic,
  94. 0,
  95. maxSize, maxSize> ComplexMatrix;
  96. const Matrix A(Matrix::Random(size, size)), B(Matrix::Random(size, size));
  97. Matrix X(size,size);
  98. const ComplexMatrix complexA(ComplexMatrix::Random(size, size));
  99. const Matrix saA = A.adjoint() * A;
  100. const Vector b(Vector::Random(size));
  101. Vector x(size);
  102. // Cholesky module
  103. Eigen::LLT<Matrix> LLT; LLT.compute(A);
  104. X = LLT.solve(B);
  105. x = LLT.solve(b);
  106. Eigen::LDLT<Matrix> LDLT; LDLT.compute(A);
  107. X = LDLT.solve(B);
  108. x = LDLT.solve(b);
  109. // Eigenvalues module
  110. Eigen::HessenbergDecomposition<ComplexMatrix> hessDecomp; hessDecomp.compute(complexA);
  111. Eigen::ComplexSchur<ComplexMatrix> cSchur(size); cSchur.compute(complexA);
  112. Eigen::ComplexEigenSolver<ComplexMatrix> cEigSolver; cEigSolver.compute(complexA);
  113. Eigen::EigenSolver<Matrix> eigSolver; eigSolver.compute(A);
  114. Eigen::SelfAdjointEigenSolver<Matrix> saEigSolver(size); saEigSolver.compute(saA);
  115. Eigen::Tridiagonalization<Matrix> tridiag; tridiag.compute(saA);
  116. // LU module
  117. Eigen::PartialPivLU<Matrix> ppLU; ppLU.compute(A);
  118. X = ppLU.solve(B);
  119. x = ppLU.solve(b);
  120. Eigen::FullPivLU<Matrix> fpLU; fpLU.compute(A);
  121. X = fpLU.solve(B);
  122. x = fpLU.solve(b);
  123. // QR module
  124. Eigen::HouseholderQR<Matrix> hQR; hQR.compute(A);
  125. X = hQR.solve(B);
  126. x = hQR.solve(b);
  127. Eigen::ColPivHouseholderQR<Matrix> cpQR; cpQR.compute(A);
  128. X = cpQR.solve(B);
  129. x = cpQR.solve(b);
  130. Eigen::FullPivHouseholderQR<Matrix> fpQR; fpQR.compute(A);
  131. // FIXME X = fpQR.solve(B);
  132. x = fpQR.solve(b);
  133. // SVD module
  134. Eigen::JacobiSVD<Matrix> jSVD; jSVD.compute(A, ComputeFullU | ComputeFullV);
  135. }
  136. void test_nomalloc()
  137. {
  138. // check that our operator new is indeed called:
  139. VERIFY_RAISES_ASSERT(MatrixXd dummy(MatrixXd::Random(3,3)));
  140. CALL_SUBTEST_1(nomalloc(Matrix<float, 1, 1>()) );
  141. CALL_SUBTEST_2(nomalloc(Matrix4d()) );
  142. CALL_SUBTEST_3(nomalloc(Matrix<float,32,32>()) );
  143. // Check decomposition modules with dynamic matrices that have a known compile-time max size (ctms)
  144. CALL_SUBTEST_4(ctms_decompositions<float>());
  145. }