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.

142 lines
5.7 KiB

  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
  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/QR>
  11. template<typename Derived1, typename Derived2>
  12. bool areNotApprox(const MatrixBase<Derived1>& m1, const MatrixBase<Derived2>& m2, typename Derived1::RealScalar epsilon = NumTraits<typename Derived1::RealScalar>::dummy_precision())
  13. {
  14. return !((m1-m2).cwiseAbs2().maxCoeff() < epsilon * epsilon
  15. * (std::max)(m1.cwiseAbs2().maxCoeff(), m2.cwiseAbs2().maxCoeff()));
  16. }
  17. template<typename MatrixType> void product(const MatrixType& m)
  18. {
  19. /* this test covers the following files:
  20. Identity.h Product.h
  21. */
  22. typedef typename MatrixType::Index Index;
  23. typedef typename MatrixType::Scalar Scalar;
  24. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> RowVectorType;
  25. typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> ColVectorType;
  26. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> RowSquareMatrixType;
  27. typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> ColSquareMatrixType;
  28. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime,
  29. MatrixType::Flags&RowMajorBit?ColMajor:RowMajor> OtherMajorMatrixType;
  30. Index rows = m.rows();
  31. Index cols = m.cols();
  32. // this test relies a lot on Random.h, and there's not much more that we can do
  33. // to test it, hence I consider that we will have tested Random.h
  34. MatrixType m1 = MatrixType::Random(rows, cols),
  35. m2 = MatrixType::Random(rows, cols),
  36. m3(rows, cols);
  37. RowSquareMatrixType
  38. identity = RowSquareMatrixType::Identity(rows, rows),
  39. square = RowSquareMatrixType::Random(rows, rows),
  40. res = RowSquareMatrixType::Random(rows, rows);
  41. ColSquareMatrixType
  42. square2 = ColSquareMatrixType::Random(cols, cols),
  43. res2 = ColSquareMatrixType::Random(cols, cols);
  44. RowVectorType v1 = RowVectorType::Random(rows);
  45. ColVectorType vc2 = ColVectorType::Random(cols), vcres(cols);
  46. OtherMajorMatrixType tm1 = m1;
  47. Scalar s1 = internal::random<Scalar>();
  48. Index r = internal::random<Index>(0, rows-1),
  49. c = internal::random<Index>(0, cols-1),
  50. c2 = internal::random<Index>(0, cols-1);
  51. // begin testing Product.h: only associativity for now
  52. // (we use Transpose.h but this doesn't count as a test for it)
  53. VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
  54. m3 = m1;
  55. m3 *= m1.transpose() * m2;
  56. VERIFY_IS_APPROX(m3, m1 * (m1.transpose()*m2));
  57. VERIFY_IS_APPROX(m3, m1 * (m1.transpose()*m2));
  58. // continue testing Product.h: distributivity
  59. VERIFY_IS_APPROX(square*(m1 + m2), square*m1+square*m2);
  60. VERIFY_IS_APPROX(square*(m1 - m2), square*m1-square*m2);
  61. // continue testing Product.h: compatibility with ScalarMultiple.h
  62. VERIFY_IS_APPROX(s1*(square*m1), (s1*square)*m1);
  63. VERIFY_IS_APPROX(s1*(square*m1), square*(m1*s1));
  64. // test Product.h together with Identity.h
  65. VERIFY_IS_APPROX(v1, identity*v1);
  66. VERIFY_IS_APPROX(v1.transpose(), v1.transpose() * identity);
  67. // again, test operator() to check const-qualification
  68. VERIFY_IS_APPROX(MatrixType::Identity(rows, cols)(r,c), static_cast<Scalar>(r==c));
  69. if (rows!=cols)
  70. VERIFY_RAISES_ASSERT(m3 = m1*m1);
  71. // test the previous tests were not screwed up because operator* returns 0
  72. // (we use the more accurate default epsilon)
  73. if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
  74. {
  75. VERIFY(areNotApprox(m1.transpose()*m2,m2.transpose()*m1));
  76. }
  77. // test optimized operator+= path
  78. res = square;
  79. res.noalias() += m1 * m2.transpose();
  80. VERIFY_IS_APPROX(res, square + m1 * m2.transpose());
  81. if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
  82. {
  83. VERIFY(areNotApprox(res,square + m2 * m1.transpose()));
  84. }
  85. vcres = vc2;
  86. vcres.noalias() += m1.transpose() * v1;
  87. VERIFY_IS_APPROX(vcres, vc2 + m1.transpose() * v1);
  88. // test optimized operator-= path
  89. res = square;
  90. res.noalias() -= m1 * m2.transpose();
  91. VERIFY_IS_APPROX(res, square - (m1 * m2.transpose()));
  92. if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
  93. {
  94. VERIFY(areNotApprox(res,square - m2 * m1.transpose()));
  95. }
  96. vcres = vc2;
  97. vcres.noalias() -= m1.transpose() * v1;
  98. VERIFY_IS_APPROX(vcres, vc2 - m1.transpose() * v1);
  99. tm1 = m1;
  100. VERIFY_IS_APPROX(tm1.transpose() * v1, m1.transpose() * v1);
  101. VERIFY_IS_APPROX(v1.transpose() * tm1, v1.transpose() * m1);
  102. // test submatrix and matrix/vector product
  103. for (int i=0; i<rows; ++i)
  104. res.row(i) = m1.row(i) * m2.transpose();
  105. VERIFY_IS_APPROX(res, m1 * m2.transpose());
  106. // the other way round:
  107. for (int i=0; i<rows; ++i)
  108. res.col(i) = m1 * m2.transpose().col(i);
  109. VERIFY_IS_APPROX(res, m1 * m2.transpose());
  110. res2 = square2;
  111. res2.noalias() += m1.transpose() * m2;
  112. VERIFY_IS_APPROX(res2, square2 + m1.transpose() * m2);
  113. if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
  114. {
  115. VERIFY(areNotApprox(res2,square2 + m2.transpose() * m1));
  116. }
  117. VERIFY_IS_APPROX(res.col(r).noalias() = square.adjoint() * square.col(r), (square.adjoint() * square.col(r)).eval());
  118. VERIFY_IS_APPROX(res.col(r).noalias() = square * square.col(r), (square * square.col(r)).eval());
  119. // inner product
  120. Scalar x = square2.row(c) * square2.col(c2);
  121. VERIFY_IS_APPROX(x, square2.row(c).transpose().cwiseProduct(square2.col(c2)).sum());
  122. }