inplace_decomposition.cpp 3.8 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2016 Gael Guennebaud <gael.guennebaud@inria.fr>
  5. //
  6. // This Source Code Form is subject to the terms of the Mozilla
  7. // Public License v. 2.0. If a copy of the MPL was not distributed
  8. // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
  9. #include "main.h"
  10. #include <Eigen/LU>
  11. #include <Eigen/Cholesky>
  12. #include <Eigen/QR>
  13. // This file test inplace decomposition through Ref<>, as supported by Cholesky, LU, and QR decompositions.
  14. template<typename DecType,typename MatrixType> void inplace(bool square = false, bool SPD = false)
  15. {
  16. typedef typename MatrixType::Scalar Scalar;
  17. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> RhsType;
  18. typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> ResType;
  19. Index rows = MatrixType::RowsAtCompileTime==Dynamic ? internal::random<Index>(2,EIGEN_TEST_MAX_SIZE/2) : Index(MatrixType::RowsAtCompileTime);
  20. Index cols = MatrixType::ColsAtCompileTime==Dynamic ? (square?rows:internal::random<Index>(2,rows)) : Index(MatrixType::ColsAtCompileTime);
  21. MatrixType A = MatrixType::Random(rows,cols);
  22. RhsType b = RhsType::Random(rows);
  23. ResType x(cols);
  24. if(SPD)
  25. {
  26. assert(square);
  27. A.topRows(cols) = A.topRows(cols).adjoint() * A.topRows(cols);
  28. A.diagonal().array() += 1e-3;
  29. }
  30. MatrixType A0 = A;
  31. MatrixType A1 = A;
  32. DecType dec(A);
  33. // Check that the content of A has been modified
  34. VERIFY_IS_NOT_APPROX( A, A0 );
  35. // Check that the decomposition is correct:
  36. if(rows==cols)
  37. {
  38. VERIFY_IS_APPROX( A0 * (x = dec.solve(b)), b );
  39. }
  40. else
  41. {
  42. VERIFY_IS_APPROX( A0.transpose() * A0 * (x = dec.solve(b)), A0.transpose() * b );
  43. }
  44. // Check that modifying A breaks the current dec:
  45. A.setRandom();
  46. if(rows==cols)
  47. {
  48. VERIFY_IS_NOT_APPROX( A0 * (x = dec.solve(b)), b );
  49. }
  50. else
  51. {
  52. VERIFY_IS_NOT_APPROX( A0.transpose() * A0 * (x = dec.solve(b)), A0.transpose() * b );
  53. }
  54. // Check that calling compute(A1) does not modify A1:
  55. A = A0;
  56. dec.compute(A1);
  57. VERIFY_IS_EQUAL(A0,A1);
  58. VERIFY_IS_NOT_APPROX( A, A0 );
  59. if(rows==cols)
  60. {
  61. VERIFY_IS_APPROX( A0 * (x = dec.solve(b)), b );
  62. }
  63. else
  64. {
  65. VERIFY_IS_APPROX( A0.transpose() * A0 * (x = dec.solve(b)), A0.transpose() * b );
  66. }
  67. }
  68. void test_inplace_decomposition()
  69. {
  70. EIGEN_UNUSED typedef Matrix<double,4,3> Matrix43d;
  71. for(int i = 0; i < g_repeat; i++) {
  72. CALL_SUBTEST_1(( inplace<LLT<Ref<MatrixXd> >, MatrixXd>(true,true) ));
  73. CALL_SUBTEST_1(( inplace<LLT<Ref<Matrix4d> >, Matrix4d>(true,true) ));
  74. CALL_SUBTEST_2(( inplace<LDLT<Ref<MatrixXd> >, MatrixXd>(true,true) ));
  75. CALL_SUBTEST_2(( inplace<LDLT<Ref<Matrix4d> >, Matrix4d>(true,true) ));
  76. CALL_SUBTEST_3(( inplace<PartialPivLU<Ref<MatrixXd> >, MatrixXd>(true,false) ));
  77. CALL_SUBTEST_3(( inplace<PartialPivLU<Ref<Matrix4d> >, Matrix4d>(true,false) ));
  78. CALL_SUBTEST_4(( inplace<FullPivLU<Ref<MatrixXd> >, MatrixXd>(true,false) ));
  79. CALL_SUBTEST_4(( inplace<FullPivLU<Ref<Matrix4d> >, Matrix4d>(true,false) ));
  80. CALL_SUBTEST_5(( inplace<HouseholderQR<Ref<MatrixXd> >, MatrixXd>(false,false) ));
  81. CALL_SUBTEST_5(( inplace<HouseholderQR<Ref<Matrix43d> >, Matrix43d>(false,false) ));
  82. CALL_SUBTEST_6(( inplace<ColPivHouseholderQR<Ref<MatrixXd> >, MatrixXd>(false,false) ));
  83. CALL_SUBTEST_6(( inplace<ColPivHouseholderQR<Ref<Matrix43d> >, Matrix43d>(false,false) ));
  84. CALL_SUBTEST_7(( inplace<FullPivHouseholderQR<Ref<MatrixXd> >, MatrixXd>(false,false) ));
  85. CALL_SUBTEST_7(( inplace<FullPivHouseholderQR<Ref<Matrix43d> >, Matrix43d>(false,false) ));
  86. CALL_SUBTEST_8(( inplace<CompleteOrthogonalDecomposition<Ref<MatrixXd> >, MatrixXd>(false,false) ));
  87. CALL_SUBTEST_8(( inplace<CompleteOrthogonalDecomposition<Ref<Matrix43d> >, Matrix43d>(false,false) ));
  88. }
  89. }