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.

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