product_trsolve.cpp 6.0 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127
  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. // destination with a non-default inner-stride
  60. // see bug 1741
  61. {
  62. typedef Matrix<Scalar,Dynamic,Dynamic> MatrixX;
  63. MatrixX buffer(2*cmRhs.rows(),2*cmRhs.cols());
  64. Map<Matrix<Scalar,Size,Cols,colmajor>,0,Stride<Dynamic,2> > map1(buffer.data(),cmRhs.rows(),cmRhs.cols(),Stride<Dynamic,2>(2*cmRhs.outerStride(),2));
  65. Map<Matrix<Scalar,Size,Cols,rowmajor>,0,Stride<Dynamic,2> > map2(buffer.data(),rmRhs.rows(),rmRhs.cols(),Stride<Dynamic,2>(2*rmRhs.outerStride(),2));
  66. buffer.setZero();
  67. VERIFY_TRSM(cmLhs.conjugate().template triangularView<Lower>(), map1);
  68. buffer.setZero();
  69. VERIFY_TRSM(cmLhs .template triangularView<Lower>(), map2);
  70. }
  71. if(Size==Dynamic)
  72. {
  73. cmLhs.resize(0,0);
  74. cmRhs.resize(0,cmRhs.cols());
  75. Matrix<Scalar,Size,Cols,colmajor> res = cmLhs.template triangularView<Lower>().solve(cmRhs);
  76. VERIFY_IS_EQUAL(res.rows(),0);
  77. VERIFY_IS_EQUAL(res.cols(),cmRhs.cols());
  78. res = cmRhs;
  79. cmLhs.template triangularView<Lower>().solveInPlace(res);
  80. VERIFY_IS_EQUAL(res.rows(),0);
  81. VERIFY_IS_EQUAL(res.cols(),cmRhs.cols());
  82. }
  83. }
  84. void test_product_trsolve()
  85. {
  86. for(int i = 0; i < g_repeat ; i++)
  87. {
  88. // matrices
  89. CALL_SUBTEST_1((trsolve<float,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  90. CALL_SUBTEST_2((trsolve<double,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  91. 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))));
  92. 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))));
  93. // vectors
  94. CALL_SUBTEST_5((trsolve<float,Dynamic,1>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  95. CALL_SUBTEST_6((trsolve<double,Dynamic,1>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  96. CALL_SUBTEST_7((trsolve<std::complex<float>,Dynamic,1>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  97. CALL_SUBTEST_8((trsolve<std::complex<double>,Dynamic,1>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))));
  98. // meta-unrollers
  99. CALL_SUBTEST_9((trsolve<float,4,1>()));
  100. CALL_SUBTEST_10((trsolve<double,4,1>()));
  101. CALL_SUBTEST_11((trsolve<std::complex<float>,4,1>()));
  102. CALL_SUBTEST_12((trsolve<float,1,1>()));
  103. CALL_SUBTEST_13((trsolve<float,1,2>()));
  104. CALL_SUBTEST_14((trsolve<float,3,1>()));
  105. }
  106. }