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.

303 lines
12 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. template<typename ArrayType> void array(const ArrayType& m)
  11. {
  12. typedef typename ArrayType::Index Index;
  13. typedef typename ArrayType::Scalar Scalar;
  14. typedef typename NumTraits<Scalar>::Real RealScalar;
  15. typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
  16. typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
  17. Index rows = m.rows();
  18. Index cols = m.cols();
  19. ArrayType m1 = ArrayType::Random(rows, cols),
  20. m2 = ArrayType::Random(rows, cols),
  21. m3(rows, cols);
  22. ColVectorType cv1 = ColVectorType::Random(rows);
  23. RowVectorType rv1 = RowVectorType::Random(cols);
  24. Scalar s1 = internal::random<Scalar>(),
  25. s2 = internal::random<Scalar>();
  26. // scalar addition
  27. VERIFY_IS_APPROX(m1 + s1, s1 + m1);
  28. VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
  29. VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
  30. VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
  31. VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
  32. VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
  33. m3 = m1;
  34. m3 += s2;
  35. VERIFY_IS_APPROX(m3, m1 + s2);
  36. m3 = m1;
  37. m3 -= s1;
  38. VERIFY_IS_APPROX(m3, m1 - s1);
  39. // scalar operators via Maps
  40. m3 = m1;
  41. ArrayType::Map(m1.data(), m1.rows(), m1.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
  42. VERIFY_IS_APPROX(m1, m3 - m2);
  43. m3 = m1;
  44. ArrayType::Map(m1.data(), m1.rows(), m1.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
  45. VERIFY_IS_APPROX(m1, m3 + m2);
  46. m3 = m1;
  47. ArrayType::Map(m1.data(), m1.rows(), m1.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
  48. VERIFY_IS_APPROX(m1, m3 * m2);
  49. m3 = m1;
  50. m2 = ArrayType::Random(rows,cols);
  51. m2 = (m2==0).select(1,m2);
  52. ArrayType::Map(m1.data(), m1.rows(), m1.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
  53. VERIFY_IS_APPROX(m1, m3 / m2);
  54. // reductions
  55. VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum());
  56. VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum());
  57. if (!internal::isApprox(m1.sum(), (m1+m2).sum(), test_precision<Scalar>()))
  58. VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
  59. VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar>()));
  60. // vector-wise ops
  61. m3 = m1;
  62. VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
  63. m3 = m1;
  64. VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
  65. m3 = m1;
  66. VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
  67. m3 = m1;
  68. VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
  69. }
  70. template<typename ArrayType> void comparisons(const ArrayType& m)
  71. {
  72. typedef typename ArrayType::Index Index;
  73. typedef typename ArrayType::Scalar Scalar;
  74. typedef typename NumTraits<Scalar>::Real RealScalar;
  75. typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> VectorType;
  76. Index rows = m.rows();
  77. Index cols = m.cols();
  78. Index r = internal::random<Index>(0, rows-1),
  79. c = internal::random<Index>(0, cols-1);
  80. ArrayType m1 = ArrayType::Random(rows, cols),
  81. m2 = ArrayType::Random(rows, cols),
  82. m3(rows, cols);
  83. VERIFY(((m1 + Scalar(1)) > m1).all());
  84. VERIFY(((m1 - Scalar(1)) < m1).all());
  85. if (rows*cols>1)
  86. {
  87. m3 = m1;
  88. m3(r,c) += 1;
  89. VERIFY(! (m1 < m3).all() );
  90. VERIFY(! (m1 > m3).all() );
  91. }
  92. // comparisons to scalar
  93. VERIFY( (m1 != (m1(r,c)+1) ).any() );
  94. VERIFY( (m1 > (m1(r,c)-1) ).any() );
  95. VERIFY( (m1 < (m1(r,c)+1) ).any() );
  96. VERIFY( (m1 == m1(r,c) ).any() );
  97. // test Select
  98. VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
  99. VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
  100. Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
  101. for (int j=0; j<cols; ++j)
  102. for (int i=0; i<rows; ++i)
  103. m3(i,j) = internal::abs(m1(i,j))<mid ? 0 : m1(i,j);
  104. VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
  105. .select(ArrayType::Zero(rows,cols),m1), m3);
  106. // shorter versions:
  107. VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
  108. .select(0,m1), m3);
  109. VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
  110. .select(m1,0), m3);
  111. // even shorter version:
  112. VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
  113. // count
  114. VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
  115. // and/or
  116. VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0);
  117. VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols);
  118. RealScalar a = m1.abs().mean();
  119. VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count());
  120. typedef Array<typename ArrayType::Index, Dynamic, 1> ArrayOfIndices;
  121. // TODO allows colwise/rowwise for array
  122. VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
  123. VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
  124. }
  125. template<typename ArrayType> void array_real(const ArrayType& m)
  126. {
  127. typedef typename ArrayType::Index Index;
  128. typedef typename ArrayType::Scalar Scalar;
  129. typedef typename NumTraits<Scalar>::Real RealScalar;
  130. Index rows = m.rows();
  131. Index cols = m.cols();
  132. ArrayType m1 = ArrayType::Random(rows, cols),
  133. m2 = ArrayType::Random(rows, cols),
  134. m3(rows, cols);
  135. Scalar s1 = internal::random<Scalar>();
  136. // these tests are mostly to check possible compilation issues.
  137. VERIFY_IS_APPROX(m1.sin(), std::sin(m1));
  138. VERIFY_IS_APPROX(m1.sin(), internal::sin(m1));
  139. VERIFY_IS_APPROX(m1.cos(), std::cos(m1));
  140. VERIFY_IS_APPROX(m1.cos(), internal::cos(m1));
  141. VERIFY_IS_APPROX(m1.asin(), std::asin(m1));
  142. VERIFY_IS_APPROX(m1.asin(), internal::asin(m1));
  143. VERIFY_IS_APPROX(m1.acos(), std::acos(m1));
  144. VERIFY_IS_APPROX(m1.acos(), internal::acos(m1));
  145. VERIFY_IS_APPROX(m1.tan(), std::tan(m1));
  146. VERIFY_IS_APPROX(m1.tan(), internal::tan(m1));
  147. VERIFY_IS_APPROX(internal::cos(m1+RealScalar(3)*m2), internal::cos((m1+RealScalar(3)*m2).eval()));
  148. VERIFY_IS_APPROX(std::cos(m1+RealScalar(3)*m2), std::cos((m1+RealScalar(3)*m2).eval()));
  149. VERIFY_IS_APPROX(m1.abs().sqrt(), std::sqrt(std::abs(m1)));
  150. VERIFY_IS_APPROX(m1.abs().sqrt(), internal::sqrt(internal::abs(m1)));
  151. VERIFY_IS_APPROX(m1.abs(), internal::sqrt(internal::abs2(m1)));
  152. VERIFY_IS_APPROX(internal::abs2(internal::real(m1)) + internal::abs2(internal::imag(m1)), internal::abs2(m1));
  153. VERIFY_IS_APPROX(internal::abs2(std::real(m1)) + internal::abs2(std::imag(m1)), internal::abs2(m1));
  154. if(!NumTraits<Scalar>::IsComplex)
  155. VERIFY_IS_APPROX(internal::real(m1), m1);
  156. VERIFY_IS_APPROX(m1.abs().log(), std::log(std::abs(m1)));
  157. VERIFY_IS_APPROX(m1.abs().log(), internal::log(internal::abs(m1)));
  158. VERIFY_IS_APPROX(m1.exp(), std::exp(m1));
  159. VERIFY_IS_APPROX(m1.exp() * m2.exp(), std::exp(m1+m2));
  160. VERIFY_IS_APPROX(m1.exp(), internal::exp(m1));
  161. VERIFY_IS_APPROX(m1.exp() / m2.exp(), std::exp(m1-m2));
  162. VERIFY_IS_APPROX(m1.pow(2), m1.square());
  163. VERIFY_IS_APPROX(std::pow(m1,2), m1.square());
  164. ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2));
  165. VERIFY_IS_APPROX(std::pow(m1,exponents), m1.square());
  166. m3 = m1.abs();
  167. VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
  168. VERIFY_IS_APPROX(std::pow(m3,RealScalar(0.5)), m3.sqrt());
  169. // scalar by array division
  170. const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon());
  171. s1 += Scalar(tiny);
  172. m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
  173. VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
  174. }
  175. template<typename ArrayType> void array_complex(const ArrayType& m)
  176. {
  177. typedef typename ArrayType::Index Index;
  178. Index rows = m.rows();
  179. Index cols = m.cols();
  180. ArrayType m1 = ArrayType::Random(rows, cols),
  181. m2(rows, cols);
  182. for (Index i = 0; i < m.rows(); ++i)
  183. for (Index j = 0; j < m.cols(); ++j)
  184. m2(i,j) = std::sqrt(m1(i,j));
  185. VERIFY_IS_APPROX(m1.sqrt(), m2);
  186. VERIFY_IS_APPROX(m1.sqrt(), std::sqrt(m1));
  187. VERIFY_IS_APPROX(m1.sqrt(), internal::sqrt(m1));
  188. }
  189. template<typename ArrayType> void min_max(const ArrayType& m)
  190. {
  191. typedef typename ArrayType::Index Index;
  192. typedef typename ArrayType::Scalar Scalar;
  193. Index rows = m.rows();
  194. Index cols = m.cols();
  195. ArrayType m1 = ArrayType::Random(rows, cols);
  196. // min/max with array
  197. Scalar maxM1 = m1.maxCoeff();
  198. Scalar minM1 = m1.minCoeff();
  199. VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1)));
  200. VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1)));
  201. VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1)));
  202. VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1)));
  203. // min/max with scalar input
  204. VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1));
  205. VERIFY_IS_APPROX(m1, (m1.min)( maxM1));
  206. VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1));
  207. VERIFY_IS_APPROX(m1, (m1.max)( minM1));
  208. }
  209. void test_array()
  210. {
  211. for(int i = 0; i < g_repeat; i++) {
  212. CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
  213. CALL_SUBTEST_2( array(Array22f()) );
  214. CALL_SUBTEST_3( array(Array44d()) );
  215. CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  216. CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  217. CALL_SUBTEST_6( array(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  218. }
  219. for(int i = 0; i < g_repeat; i++) {
  220. CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
  221. CALL_SUBTEST_2( comparisons(Array22f()) );
  222. CALL_SUBTEST_3( comparisons(Array44d()) );
  223. CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  224. CALL_SUBTEST_6( comparisons(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  225. }
  226. for(int i = 0; i < g_repeat; i++) {
  227. CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) );
  228. CALL_SUBTEST_2( min_max(Array22f()) );
  229. CALL_SUBTEST_3( min_max(Array44d()) );
  230. CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  231. CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  232. }
  233. for(int i = 0; i < g_repeat; i++) {
  234. CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
  235. CALL_SUBTEST_2( array_real(Array22f()) );
  236. CALL_SUBTEST_3( array_real(Array44d()) );
  237. CALL_SUBTEST_5( array_real(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  238. }
  239. for(int i = 0; i < g_repeat; i++) {
  240. CALL_SUBTEST_4( array_complex(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  241. }
  242. VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
  243. VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
  244. VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
  245. typedef CwiseUnaryOp<internal::scalar_sum_op<double>, ArrayXd > Xpr;
  246. VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
  247. ArrayBase<Xpr>
  248. >::value));
  249. }