array_for_matrix.cpp 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300
  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. template<typename MatrixType> void array_for_matrix(const MatrixType& m)
  11. {
  12. typedef typename MatrixType::Scalar Scalar;
  13. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> ColVectorType;
  14. typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
  15. Index rows = m.rows();
  16. Index cols = m.cols();
  17. MatrixType m1 = MatrixType::Random(rows, cols),
  18. m2 = MatrixType::Random(rows, cols),
  19. m3(rows, cols);
  20. ColVectorType cv1 = ColVectorType::Random(rows);
  21. RowVectorType rv1 = RowVectorType::Random(cols);
  22. Scalar s1 = internal::random<Scalar>(),
  23. s2 = internal::random<Scalar>();
  24. // scalar addition
  25. VERIFY_IS_APPROX(m1.array() + s1, s1 + m1.array());
  26. VERIFY_IS_APPROX((m1.array() + s1).matrix(), MatrixType::Constant(rows,cols,s1) + m1);
  27. VERIFY_IS_APPROX(((m1*Scalar(2)).array() - s2).matrix(), (m1+m1) - MatrixType::Constant(rows,cols,s2) );
  28. m3 = m1;
  29. m3.array() += s2;
  30. VERIFY_IS_APPROX(m3, (m1.array() + s2).matrix());
  31. m3 = m1;
  32. m3.array() -= s1;
  33. VERIFY_IS_APPROX(m3, (m1.array() - s1).matrix());
  34. // reductions
  35. VERIFY_IS_MUCH_SMALLER_THAN(m1.colwise().sum().sum() - m1.sum(), m1.squaredNorm());
  36. VERIFY_IS_MUCH_SMALLER_THAN(m1.rowwise().sum().sum() - m1.sum(), m1.squaredNorm());
  37. VERIFY_IS_MUCH_SMALLER_THAN(m1.colwise().sum() + m2.colwise().sum() - (m1+m2).colwise().sum(), (m1+m2).squaredNorm());
  38. VERIFY_IS_MUCH_SMALLER_THAN(m1.rowwise().sum() - m2.rowwise().sum() - (m1-m2).rowwise().sum(), (m1-m2).squaredNorm());
  39. VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
  40. // vector-wise ops
  41. m3 = m1;
  42. VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
  43. m3 = m1;
  44. VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
  45. m3 = m1;
  46. VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
  47. m3 = m1;
  48. VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
  49. // empty objects
  50. VERIFY_IS_APPROX(m1.block(0,0,0,cols).colwise().sum(), RowVectorType::Zero(cols));
  51. VERIFY_IS_APPROX(m1.block(0,0,rows,0).rowwise().prod(), ColVectorType::Ones(rows));
  52. // verify the const accessors exist
  53. const Scalar& ref_m1 = m.matrix().array().coeffRef(0);
  54. const Scalar& ref_m2 = m.matrix().array().coeffRef(0,0);
  55. const Scalar& ref_a1 = m.array().matrix().coeffRef(0);
  56. const Scalar& ref_a2 = m.array().matrix().coeffRef(0,0);
  57. VERIFY(&ref_a1 == &ref_m1);
  58. VERIFY(&ref_a2 == &ref_m2);
  59. // Check write accessors:
  60. m1.array().coeffRef(0,0) = 1;
  61. VERIFY_IS_APPROX(m1(0,0),Scalar(1));
  62. m1.array()(0,0) = 2;
  63. VERIFY_IS_APPROX(m1(0,0),Scalar(2));
  64. m1.array().matrix().coeffRef(0,0) = 3;
  65. VERIFY_IS_APPROX(m1(0,0),Scalar(3));
  66. m1.array().matrix()(0,0) = 4;
  67. VERIFY_IS_APPROX(m1(0,0),Scalar(4));
  68. }
  69. template<typename MatrixType> void comparisons(const MatrixType& m)
  70. {
  71. using std::abs;
  72. typedef typename MatrixType::Scalar Scalar;
  73. typedef typename NumTraits<Scalar>::Real RealScalar;
  74. Index rows = m.rows();
  75. Index cols = m.cols();
  76. Index r = internal::random<Index>(0, rows-1),
  77. c = internal::random<Index>(0, cols-1);
  78. MatrixType m1 = MatrixType::Random(rows, cols),
  79. m2 = MatrixType::Random(rows, cols),
  80. m3(rows, cols);
  81. VERIFY(((m1.array() + Scalar(1)) > m1.array()).all());
  82. VERIFY(((m1.array() - Scalar(1)) < m1.array()).all());
  83. if (rows*cols>1)
  84. {
  85. m3 = m1;
  86. m3(r,c) += 1;
  87. VERIFY(! (m1.array() < m3.array()).all() );
  88. VERIFY(! (m1.array() > m3.array()).all() );
  89. }
  90. // comparisons to scalar
  91. VERIFY( (m1.array() != (m1(r,c)+1) ).any() );
  92. VERIFY( (m1.array() > (m1(r,c)-1) ).any() );
  93. VERIFY( (m1.array() < (m1(r,c)+1) ).any() );
  94. VERIFY( (m1.array() == m1(r,c) ).any() );
  95. VERIFY( m1.cwiseEqual(m1(r,c)).any() );
  96. // test Select
  97. VERIFY_IS_APPROX( (m1.array()<m2.array()).select(m1,m2), m1.cwiseMin(m2) );
  98. VERIFY_IS_APPROX( (m1.array()>m2.array()).select(m1,m2), m1.cwiseMax(m2) );
  99. Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
  100. for (int j=0; j<cols; ++j)
  101. for (int i=0; i<rows; ++i)
  102. m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
  103. VERIFY_IS_APPROX( (m1.array().abs()<MatrixType::Constant(rows,cols,mid).array())
  104. .select(MatrixType::Zero(rows,cols),m1), m3);
  105. // shorter versions:
  106. VERIFY_IS_APPROX( (m1.array().abs()<MatrixType::Constant(rows,cols,mid).array())
  107. .select(0,m1), m3);
  108. VERIFY_IS_APPROX( (m1.array().abs()>=MatrixType::Constant(rows,cols,mid).array())
  109. .select(m1,0), m3);
  110. // even shorter version:
  111. VERIFY_IS_APPROX( (m1.array().abs()<mid).select(0,m1), m3);
  112. // count
  113. VERIFY(((m1.array().abs()+1)>RealScalar(0.1)).count() == rows*cols);
  114. // and/or
  115. VERIFY( ((m1.array()<RealScalar(0)).matrix() && (m1.array()>RealScalar(0)).matrix()).count() == 0);
  116. VERIFY( ((m1.array()<RealScalar(0)).matrix() || (m1.array()>=RealScalar(0)).matrix()).count() == rows*cols);
  117. RealScalar a = m1.cwiseAbs().mean();
  118. VERIFY( ((m1.array()<-a).matrix() || (m1.array()>a).matrix()).count() == (m1.cwiseAbs().array()>a).count());
  119. typedef Matrix<typename MatrixType::Index, Dynamic, 1> VectorOfIndices;
  120. // TODO allows colwise/rowwise for array
  121. VERIFY_IS_APPROX(((m1.array().abs()+1)>RealScalar(0.1)).matrix().colwise().count(), VectorOfIndices::Constant(cols,rows).transpose());
  122. VERIFY_IS_APPROX(((m1.array().abs()+1)>RealScalar(0.1)).matrix().rowwise().count(), VectorOfIndices::Constant(rows, cols));
  123. }
  124. template<typename VectorType> void lpNorm(const VectorType& v)
  125. {
  126. using std::sqrt;
  127. typedef typename VectorType::RealScalar RealScalar;
  128. VectorType u = VectorType::Random(v.size());
  129. if(v.size()==0)
  130. {
  131. VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), RealScalar(0));
  132. VERIFY_IS_APPROX(u.template lpNorm<1>(), RealScalar(0));
  133. VERIFY_IS_APPROX(u.template lpNorm<2>(), RealScalar(0));
  134. VERIFY_IS_APPROX(u.template lpNorm<5>(), RealScalar(0));
  135. }
  136. else
  137. {
  138. VERIFY_IS_APPROX(u.template lpNorm<Infinity>(), u.cwiseAbs().maxCoeff());
  139. }
  140. VERIFY_IS_APPROX(u.template lpNorm<1>(), u.cwiseAbs().sum());
  141. VERIFY_IS_APPROX(u.template lpNorm<2>(), sqrt(u.array().abs().square().sum()));
  142. VERIFY_IS_APPROX(numext::pow(u.template lpNorm<5>(), typename VectorType::RealScalar(5)), u.array().abs().pow(5).sum());
  143. }
  144. template<typename MatrixType> void cwise_min_max(const MatrixType& m)
  145. {
  146. typedef typename MatrixType::Scalar Scalar;
  147. Index rows = m.rows();
  148. Index cols = m.cols();
  149. MatrixType m1 = MatrixType::Random(rows, cols);
  150. // min/max with array
  151. Scalar maxM1 = m1.maxCoeff();
  152. Scalar minM1 = m1.minCoeff();
  153. VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, minM1), m1.cwiseMin(MatrixType::Constant(rows,cols, minM1)));
  154. VERIFY_IS_APPROX(m1, m1.cwiseMin(MatrixType::Constant(rows,cols, maxM1)));
  155. VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, maxM1), m1.cwiseMax(MatrixType::Constant(rows,cols, maxM1)));
  156. VERIFY_IS_APPROX(m1, m1.cwiseMax(MatrixType::Constant(rows,cols, minM1)));
  157. // min/max with scalar input
  158. VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, minM1), m1.cwiseMin( minM1));
  159. VERIFY_IS_APPROX(m1, m1.cwiseMin(maxM1));
  160. VERIFY_IS_APPROX(-m1, (-m1).cwiseMin(-minM1));
  161. VERIFY_IS_APPROX(-m1.array(), ((-m1).array().min)( -minM1));
  162. VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, maxM1), m1.cwiseMax( maxM1));
  163. VERIFY_IS_APPROX(m1, m1.cwiseMax(minM1));
  164. VERIFY_IS_APPROX(-m1, (-m1).cwiseMax(-maxM1));
  165. VERIFY_IS_APPROX(-m1.array(), ((-m1).array().max)(-maxM1));
  166. VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, minM1).array(), (m1.array().min)( minM1));
  167. VERIFY_IS_APPROX(m1.array(), (m1.array().min)( maxM1));
  168. VERIFY_IS_APPROX(MatrixType::Constant(rows,cols, maxM1).array(), (m1.array().max)( maxM1));
  169. VERIFY_IS_APPROX(m1.array(), (m1.array().max)( minM1));
  170. }
  171. template<typename MatrixTraits> void resize(const MatrixTraits& t)
  172. {
  173. typedef typename MatrixTraits::Scalar Scalar;
  174. typedef Matrix<Scalar,Dynamic,Dynamic> MatrixType;
  175. typedef Array<Scalar,Dynamic,Dynamic> Array2DType;
  176. typedef Matrix<Scalar,Dynamic,1> VectorType;
  177. typedef Array<Scalar,Dynamic,1> Array1DType;
  178. Index rows = t.rows(), cols = t.cols();
  179. MatrixType m(rows,cols);
  180. VectorType v(rows);
  181. Array2DType a2(rows,cols);
  182. Array1DType a1(rows);
  183. m.array().resize(rows+1,cols+1);
  184. VERIFY(m.rows()==rows+1 && m.cols()==cols+1);
  185. a2.matrix().resize(rows+1,cols+1);
  186. VERIFY(a2.rows()==rows+1 && a2.cols()==cols+1);
  187. v.array().resize(cols);
  188. VERIFY(v.size()==cols);
  189. a1.matrix().resize(cols);
  190. VERIFY(a1.size()==cols);
  191. }
  192. template<int>
  193. void regression_bug_654()
  194. {
  195. ArrayXf a = RowVectorXf(3);
  196. VectorXf v = Array<float,1,Dynamic>(3);
  197. }
  198. // Check propagation of LvalueBit through Array/Matrix-Wrapper
  199. template<int>
  200. void regrrssion_bug_1410()
  201. {
  202. const Matrix4i M;
  203. const Array4i A;
  204. ArrayWrapper<const Matrix4i> MA = M.array();
  205. MA.row(0);
  206. MatrixWrapper<const Array4i> AM = A.matrix();
  207. AM.row(0);
  208. VERIFY((internal::traits<ArrayWrapper<const Matrix4i> >::Flags&LvalueBit)==0);
  209. VERIFY((internal::traits<MatrixWrapper<const Array4i> >::Flags&LvalueBit)==0);
  210. VERIFY((internal::traits<ArrayWrapper<Matrix4i> >::Flags&LvalueBit)==LvalueBit);
  211. VERIFY((internal::traits<MatrixWrapper<Array4i> >::Flags&LvalueBit)==LvalueBit);
  212. }
  213. void test_array_for_matrix()
  214. {
  215. for(int i = 0; i < g_repeat; i++) {
  216. CALL_SUBTEST_1( array_for_matrix(Matrix<float, 1, 1>()) );
  217. CALL_SUBTEST_2( array_for_matrix(Matrix2f()) );
  218. CALL_SUBTEST_3( array_for_matrix(Matrix4d()) );
  219. CALL_SUBTEST_4( array_for_matrix(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  220. CALL_SUBTEST_5( array_for_matrix(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  221. CALL_SUBTEST_6( array_for_matrix(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  222. }
  223. for(int i = 0; i < g_repeat; i++) {
  224. CALL_SUBTEST_1( comparisons(Matrix<float, 1, 1>()) );
  225. CALL_SUBTEST_2( comparisons(Matrix2f()) );
  226. CALL_SUBTEST_3( comparisons(Matrix4d()) );
  227. CALL_SUBTEST_5( comparisons(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  228. CALL_SUBTEST_6( comparisons(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  229. }
  230. for(int i = 0; i < g_repeat; i++) {
  231. CALL_SUBTEST_1( cwise_min_max(Matrix<float, 1, 1>()) );
  232. CALL_SUBTEST_2( cwise_min_max(Matrix2f()) );
  233. CALL_SUBTEST_3( cwise_min_max(Matrix4d()) );
  234. CALL_SUBTEST_5( cwise_min_max(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  235. CALL_SUBTEST_6( cwise_min_max(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  236. }
  237. for(int i = 0; i < g_repeat; i++) {
  238. CALL_SUBTEST_1( lpNorm(Matrix<float, 1, 1>()) );
  239. CALL_SUBTEST_2( lpNorm(Vector2f()) );
  240. CALL_SUBTEST_7( lpNorm(Vector3d()) );
  241. CALL_SUBTEST_8( lpNorm(Vector4f()) );
  242. CALL_SUBTEST_5( lpNorm(VectorXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  243. CALL_SUBTEST_4( lpNorm(VectorXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  244. }
  245. CALL_SUBTEST_5( lpNorm(VectorXf(0)) );
  246. CALL_SUBTEST_4( lpNorm(VectorXcf(0)) );
  247. for(int i = 0; i < g_repeat; i++) {
  248. CALL_SUBTEST_4( resize(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  249. CALL_SUBTEST_5( resize(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  250. CALL_SUBTEST_6( resize(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  251. }
  252. CALL_SUBTEST_6( regression_bug_654<0>() );
  253. CALL_SUBTEST_6( regrrssion_bug_1410<0>() );
  254. }