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.

101 lines
4.9 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2008-2009 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. #define VERIFY_TRSM(TRI,XB) { \
  11. (XB).setRandom(); ref = (XB); \
  12. (TRI).solveInPlace(XB); \
  13. VERIFY_IS_APPROX((TRI).toDenseMatrix() * (XB), ref); \
  14. (XB).setRandom(); ref = (XB); \
  15. (XB) = (TRI).solve(XB); \
  16. VERIFY_IS_APPROX((TRI).toDenseMatrix() * (XB), ref); \
  17. }
  18. #define VERIFY_TRSM_ONTHERIGHT(TRI,XB) { \
  19. (XB).setRandom(); ref = (XB); \
  20. (TRI).transpose().template solveInPlace<OnTheRight>(XB.transpose()); \
  21. VERIFY_IS_APPROX((XB).transpose() * (TRI).transpose().toDenseMatrix(), ref.transpose()); \
  22. (XB).setRandom(); ref = (XB); \
  23. (XB).transpose() = (TRI).transpose().template solve<OnTheRight>(XB.transpose()); \
  24. VERIFY_IS_APPROX((XB).transpose() * (TRI).transpose().toDenseMatrix(), ref.transpose()); \
  25. }
  26. template<typename Scalar,int Size, int Cols> void trsolve(int size=Size,int cols=Cols)
  27. {
  28. typedef typename NumTraits<Scalar>::Real RealScalar;
  29. Matrix<Scalar,Size,Size,ColMajor> cmLhs(size,size);
  30. Matrix<Scalar,Size,Size,RowMajor> rmLhs(size,size);
  31. enum { colmajor = Size==1 ? RowMajor : ColMajor,
  32. rowmajor = Cols==1 ? ColMajor : RowMajor };
  33. Matrix<Scalar,Size,Cols,colmajor> cmRhs(size,cols);
  34. Matrix<Scalar,Size,Cols,rowmajor> rmRhs(size,cols);
  35. Matrix<Scalar,Dynamic,Dynamic,colmajor> ref(size,cols);
  36. cmLhs.setRandom(); cmLhs *= static_cast<RealScalar>(0.1); cmLhs.diagonal().array() += static_cast<RealScalar>(1);
  37. rmLhs.setRandom(); rmLhs *= static_cast<RealScalar>(0.1); rmLhs.diagonal().array() += static_cast<RealScalar>(1);
  38. VERIFY_TRSM(cmLhs.conjugate().template triangularView<Lower>(), cmRhs);
  39. VERIFY_TRSM(cmLhs.adjoint() .template triangularView<Lower>(), cmRhs);
  40. VERIFY_TRSM(cmLhs .template triangularView<Upper>(), cmRhs);
  41. VERIFY_TRSM(cmLhs .template triangularView<Lower>(), rmRhs);
  42. VERIFY_TRSM(cmLhs.conjugate().template triangularView<Upper>(), rmRhs);
  43. VERIFY_TRSM(cmLhs.adjoint() .template triangularView<Upper>(), rmRhs);
  44. VERIFY_TRSM(cmLhs.conjugate().template triangularView<UnitLower>(), cmRhs);
  45. VERIFY_TRSM(cmLhs .template triangularView<UnitUpper>(), rmRhs);
  46. VERIFY_TRSM(rmLhs .template triangularView<Lower>(), cmRhs);
  47. VERIFY_TRSM(rmLhs.conjugate().template triangularView<UnitUpper>(), rmRhs);
  48. VERIFY_TRSM_ONTHERIGHT(cmLhs.conjugate().template triangularView<Lower>(), cmRhs);
  49. VERIFY_TRSM_ONTHERIGHT(cmLhs .template triangularView<Upper>(), cmRhs);
  50. VERIFY_TRSM_ONTHERIGHT(cmLhs .template triangularView<Lower>(), rmRhs);
  51. VERIFY_TRSM_ONTHERIGHT(cmLhs.conjugate().template triangularView<Upper>(), rmRhs);
  52. VERIFY_TRSM_ONTHERIGHT(cmLhs.conjugate().template triangularView<UnitLower>(), cmRhs);
  53. VERIFY_TRSM_ONTHERIGHT(cmLhs .template triangularView<UnitUpper>(), rmRhs);
  54. VERIFY_TRSM_ONTHERIGHT(rmLhs .template triangularView<Lower>(), cmRhs);
  55. VERIFY_TRSM_ONTHERIGHT(rmLhs.conjugate().template triangularView<UnitUpper>(), rmRhs);
  56. int c = internal::random<int>(0,cols-1);
  57. VERIFY_TRSM(rmLhs.template triangularView<Lower>(), rmRhs.col(c));
  58. VERIFY_TRSM(cmLhs.template triangularView<Lower>(), rmRhs.col(c));
  59. }
  60. void test_product_trsolve()
  61. {
  62. for(int i = 0; i < g_repeat ; i++)
  63. {
  64. // matrices
  65. CALL_SUBTEST_1((trsolve<float,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  66. CALL_SUBTEST_2((trsolve<double,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  67. CALL_SUBTEST_3((trsolve<std::complex<float>,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2),internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))));
  68. CALL_SUBTEST_4((trsolve<std::complex<double>,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2),internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))));
  69. // vectors
  70. CALL_SUBTEST_5((trsolve<float,Dynamic,1>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  71. CALL_SUBTEST_6((trsolve<double,Dynamic,1>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  72. CALL_SUBTEST_7((trsolve<std::complex<float>,Dynamic,1>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  73. CALL_SUBTEST_8((trsolve<std::complex<double>,Dynamic,1>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  74. // meta-unrollers
  75. CALL_SUBTEST_9((trsolve<float,4,1>()));
  76. CALL_SUBTEST_10((trsolve<double,4,1>()));
  77. CALL_SUBTEST_11((trsolve<std::complex<float>,4,1>()));
  78. CALL_SUBTEST_12((trsolve<float,1,1>()));
  79. CALL_SUBTEST_13((trsolve<float,1,2>()));
  80. CALL_SUBTEST_14((trsolve<float,3,1>()));
  81. }
  82. }