product.h 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257
  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::Scalar Scalar;
  23. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> RowVectorType;
  24. typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> ColVectorType;
  25. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> RowSquareMatrixType;
  26. typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> ColSquareMatrixType;
  27. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime,
  28. MatrixType::Flags&RowMajorBit?ColMajor:RowMajor> OtherMajorMatrixType;
  29. Index rows = m.rows();
  30. Index cols = m.cols();
  31. // this test relies a lot on Random.h, and there's not much more that we can do
  32. // to test it, hence I consider that we will have tested Random.h
  33. MatrixType m1 = MatrixType::Random(rows, cols),
  34. m2 = MatrixType::Random(rows, cols),
  35. m3(rows, cols);
  36. RowSquareMatrixType
  37. identity = RowSquareMatrixType::Identity(rows, rows),
  38. square = RowSquareMatrixType::Random(rows, rows),
  39. res = RowSquareMatrixType::Random(rows, rows);
  40. ColSquareMatrixType
  41. square2 = ColSquareMatrixType::Random(cols, cols),
  42. res2 = ColSquareMatrixType::Random(cols, cols);
  43. RowVectorType v1 = RowVectorType::Random(rows);
  44. ColVectorType vc2 = ColVectorType::Random(cols), vcres(cols);
  45. OtherMajorMatrixType tm1 = m1;
  46. Scalar s1 = internal::random<Scalar>();
  47. Index r = internal::random<Index>(0, rows-1),
  48. c = internal::random<Index>(0, cols-1),
  49. c2 = internal::random<Index>(0, cols-1);
  50. // begin testing Product.h: only associativity for now
  51. // (we use Transpose.h but this doesn't count as a test for it)
  52. VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
  53. m3 = m1;
  54. m3 *= m1.transpose() * m2;
  55. VERIFY_IS_APPROX(m3, m1 * (m1.transpose()*m2));
  56. VERIFY_IS_APPROX(m3, m1 * (m1.transpose()*m2));
  57. // continue testing Product.h: distributivity
  58. VERIFY_IS_APPROX(square*(m1 + m2), square*m1+square*m2);
  59. VERIFY_IS_APPROX(square*(m1 - m2), square*m1-square*m2);
  60. // continue testing Product.h: compatibility with ScalarMultiple.h
  61. VERIFY_IS_APPROX(s1*(square*m1), (s1*square)*m1);
  62. VERIFY_IS_APPROX(s1*(square*m1), square*(m1*s1));
  63. // test Product.h together with Identity.h
  64. VERIFY_IS_APPROX(v1, identity*v1);
  65. VERIFY_IS_APPROX(v1.transpose(), v1.transpose() * identity);
  66. // again, test operator() to check const-qualification
  67. VERIFY_IS_APPROX(MatrixType::Identity(rows, cols)(r,c), static_cast<Scalar>(r==c));
  68. if (rows!=cols)
  69. VERIFY_RAISES_ASSERT(m3 = m1*m1);
  70. // test the previous tests were not screwed up because operator* returns 0
  71. // (we use the more accurate default epsilon)
  72. if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
  73. {
  74. VERIFY(areNotApprox(m1.transpose()*m2,m2.transpose()*m1));
  75. }
  76. // test optimized operator+= path
  77. res = square;
  78. res.noalias() += m1 * m2.transpose();
  79. VERIFY_IS_APPROX(res, square + m1 * m2.transpose());
  80. if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
  81. {
  82. VERIFY(areNotApprox(res,square + m2 * m1.transpose()));
  83. }
  84. vcres = vc2;
  85. vcres.noalias() += m1.transpose() * v1;
  86. VERIFY_IS_APPROX(vcres, vc2 + m1.transpose() * v1);
  87. // test optimized operator-= path
  88. res = square;
  89. res.noalias() -= m1 * m2.transpose();
  90. VERIFY_IS_APPROX(res, square - (m1 * m2.transpose()));
  91. if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
  92. {
  93. VERIFY(areNotApprox(res,square - m2 * m1.transpose()));
  94. }
  95. vcres = vc2;
  96. vcres.noalias() -= m1.transpose() * v1;
  97. VERIFY_IS_APPROX(vcres, vc2 - m1.transpose() * v1);
  98. // test scaled products
  99. res = square;
  100. res.noalias() = s1 * m1 * m2.transpose();
  101. VERIFY_IS_APPROX(res, ((s1*m1).eval() * m2.transpose()));
  102. res = square;
  103. res.noalias() += s1 * m1 * m2.transpose();
  104. VERIFY_IS_APPROX(res, square + ((s1*m1).eval() * m2.transpose()));
  105. res = square;
  106. res.noalias() -= s1 * m1 * m2.transpose();
  107. VERIFY_IS_APPROX(res, square - ((s1*m1).eval() * m2.transpose()));
  108. // test d ?= a+b*c rules
  109. res.noalias() = square + m1 * m2.transpose();
  110. VERIFY_IS_APPROX(res, square + m1 * m2.transpose());
  111. res.noalias() += square + m1 * m2.transpose();
  112. VERIFY_IS_APPROX(res, 2*(square + m1 * m2.transpose()));
  113. res.noalias() -= square + m1 * m2.transpose();
  114. VERIFY_IS_APPROX(res, square + m1 * m2.transpose());
  115. // test d ?= a-b*c rules
  116. res.noalias() = square - m1 * m2.transpose();
  117. VERIFY_IS_APPROX(res, square - m1 * m2.transpose());
  118. res.noalias() += square - m1 * m2.transpose();
  119. VERIFY_IS_APPROX(res, 2*(square - m1 * m2.transpose()));
  120. res.noalias() -= square - m1 * m2.transpose();
  121. VERIFY_IS_APPROX(res, square - m1 * m2.transpose());
  122. tm1 = m1;
  123. VERIFY_IS_APPROX(tm1.transpose() * v1, m1.transpose() * v1);
  124. VERIFY_IS_APPROX(v1.transpose() * tm1, v1.transpose() * m1);
  125. // test submatrix and matrix/vector product
  126. for (int i=0; i<rows; ++i)
  127. res.row(i) = m1.row(i) * m2.transpose();
  128. VERIFY_IS_APPROX(res, m1 * m2.transpose());
  129. // the other way round:
  130. for (int i=0; i<rows; ++i)
  131. res.col(i) = m1 * m2.transpose().col(i);
  132. VERIFY_IS_APPROX(res, m1 * m2.transpose());
  133. res2 = square2;
  134. res2.noalias() += m1.transpose() * m2;
  135. VERIFY_IS_APPROX(res2, square2 + m1.transpose() * m2);
  136. if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
  137. {
  138. VERIFY(areNotApprox(res2,square2 + m2.transpose() * m1));
  139. }
  140. VERIFY_IS_APPROX(res.col(r).noalias() = square.adjoint() * square.col(r), (square.adjoint() * square.col(r)).eval());
  141. VERIFY_IS_APPROX(res.col(r).noalias() = square * square.col(r), (square * square.col(r)).eval());
  142. // vector at runtime (see bug 1166)
  143. {
  144. RowSquareMatrixType ref(square);
  145. ColSquareMatrixType ref2(square2);
  146. ref = res = square;
  147. VERIFY_IS_APPROX(res.block(0,0,1,rows).noalias() = m1.col(0).transpose() * square.transpose(), (ref.row(0) = m1.col(0).transpose() * square.transpose()));
  148. VERIFY_IS_APPROX(res.block(0,0,1,rows).noalias() = m1.block(0,0,rows,1).transpose() * square.transpose(), (ref.row(0) = m1.col(0).transpose() * square.transpose()));
  149. VERIFY_IS_APPROX(res.block(0,0,1,rows).noalias() = m1.col(0).transpose() * square, (ref.row(0) = m1.col(0).transpose() * square));
  150. VERIFY_IS_APPROX(res.block(0,0,1,rows).noalias() = m1.block(0,0,rows,1).transpose() * square, (ref.row(0) = m1.col(0).transpose() * square));
  151. ref2 = res2 = square2;
  152. VERIFY_IS_APPROX(res2.block(0,0,1,cols).noalias() = m1.row(0) * square2.transpose(), (ref2.row(0) = m1.row(0) * square2.transpose()));
  153. VERIFY_IS_APPROX(res2.block(0,0,1,cols).noalias() = m1.block(0,0,1,cols) * square2.transpose(), (ref2.row(0) = m1.row(0) * square2.transpose()));
  154. VERIFY_IS_APPROX(res2.block(0,0,1,cols).noalias() = m1.row(0) * square2, (ref2.row(0) = m1.row(0) * square2));
  155. VERIFY_IS_APPROX(res2.block(0,0,1,cols).noalias() = m1.block(0,0,1,cols) * square2, (ref2.row(0) = m1.row(0) * square2));
  156. }
  157. // vector.block() (see bug 1283)
  158. {
  159. RowVectorType w1(rows);
  160. VERIFY_IS_APPROX(square * v1.block(0,0,rows,1), square * v1);
  161. VERIFY_IS_APPROX(w1.noalias() = square * v1.block(0,0,rows,1), square * v1);
  162. VERIFY_IS_APPROX(w1.block(0,0,rows,1).noalias() = square * v1.block(0,0,rows,1), square * v1);
  163. Matrix<Scalar,1,MatrixType::ColsAtCompileTime> w2(cols);
  164. VERIFY_IS_APPROX(vc2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
  165. VERIFY_IS_APPROX(w2.noalias() = vc2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
  166. VERIFY_IS_APPROX(w2.block(0,0,1,cols).noalias() = vc2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
  167. vc2 = square2.block(0,0,1,cols).transpose();
  168. VERIFY_IS_APPROX(square2.block(0,0,1,cols) * square2, vc2.transpose() * square2);
  169. VERIFY_IS_APPROX(w2.noalias() = square2.block(0,0,1,cols) * square2, vc2.transpose() * square2);
  170. VERIFY_IS_APPROX(w2.block(0,0,1,cols).noalias() = square2.block(0,0,1,cols) * square2, vc2.transpose() * square2);
  171. vc2 = square2.block(0,0,cols,1);
  172. VERIFY_IS_APPROX(square2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
  173. VERIFY_IS_APPROX(w2.noalias() = square2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
  174. VERIFY_IS_APPROX(w2.block(0,0,1,cols).noalias() = square2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
  175. }
  176. // inner product
  177. {
  178. Scalar x = square2.row(c) * square2.col(c2);
  179. VERIFY_IS_APPROX(x, square2.row(c).transpose().cwiseProduct(square2.col(c2)).sum());
  180. }
  181. // outer product
  182. {
  183. VERIFY_IS_APPROX(m1.col(c) * m1.row(r), m1.block(0,c,rows,1) * m1.block(r,0,1,cols));
  184. VERIFY_IS_APPROX(m1.row(r).transpose() * m1.col(c).transpose(), m1.block(r,0,1,cols).transpose() * m1.block(0,c,rows,1).transpose());
  185. VERIFY_IS_APPROX(m1.block(0,c,rows,1) * m1.row(r), m1.block(0,c,rows,1) * m1.block(r,0,1,cols));
  186. VERIFY_IS_APPROX(m1.col(c) * m1.block(r,0,1,cols), m1.block(0,c,rows,1) * m1.block(r,0,1,cols));
  187. VERIFY_IS_APPROX(m1.leftCols(1) * m1.row(r), m1.block(0,0,rows,1) * m1.block(r,0,1,cols));
  188. VERIFY_IS_APPROX(m1.col(c) * m1.topRows(1), m1.block(0,c,rows,1) * m1.block(0,0,1,cols));
  189. }
  190. // Aliasing
  191. {
  192. ColVectorType x(cols); x.setRandom();
  193. ColVectorType z(x);
  194. ColVectorType y(cols); y.setZero();
  195. ColSquareMatrixType A(cols,cols); A.setRandom();
  196. // CwiseBinaryOp
  197. VERIFY_IS_APPROX(x = y + A*x, A*z);
  198. x = z;
  199. // CwiseUnaryOp
  200. VERIFY_IS_APPROX(x = Scalar(1.)*(A*x), A*z);
  201. }
  202. // regression for blas_trais
  203. {
  204. VERIFY_IS_APPROX(square * (square*square).transpose(), square * square.transpose() * square.transpose());
  205. VERIFY_IS_APPROX(square * (-(square*square)), -square * square * square);
  206. VERIFY_IS_APPROX(square * (s1*(square*square)), s1 * square * square * square);
  207. VERIFY_IS_APPROX(square * (square*square).conjugate(), square * square.conjugate() * square.conjugate());
  208. }
  209. // destination with a non-default inner-stride
  210. // see bug 1741
  211. if(!MatrixType::IsRowMajor)
  212. {
  213. typedef Matrix<Scalar,Dynamic,Dynamic> MatrixX;
  214. MatrixX buffer(2*rows,2*rows);
  215. Map<RowSquareMatrixType,0,Stride<Dynamic,2> > map1(buffer.data(),rows,rows,Stride<Dynamic,2>(2*rows,2));
  216. buffer.setZero();
  217. VERIFY_IS_APPROX(map1 = m1 * m2.transpose(), (m1 * m2.transpose()).eval());
  218. buffer.setZero();
  219. VERIFY_IS_APPROX(map1.noalias() = m1 * m2.transpose(), (m1 * m2.transpose()).eval());
  220. buffer.setZero();
  221. VERIFY_IS_APPROX(map1.noalias() += m1 * m2.transpose(), (m1 * m2.transpose()).eval());
  222. }
  223. }