basicstuff.cpp 13 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356
  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. #define EIGEN_NO_STATIC_ASSERT
  10. #include "main.h"
  11. #include "random_without_cast_overflow.h"
  12. template<typename MatrixType> void basicStuff(const MatrixType& m)
  13. {
  14. typedef typename MatrixType::Scalar Scalar;
  15. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
  16. typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
  17. Index rows = m.rows();
  18. Index cols = m.cols();
  19. // this test relies a lot on Random.h, and there's not much more that we can do
  20. // to test it, hence I consider that we will have tested Random.h
  21. MatrixType m1 = MatrixType::Random(rows, cols),
  22. m2 = MatrixType::Random(rows, cols),
  23. m3(rows, cols),
  24. mzero = MatrixType::Zero(rows, cols),
  25. square = Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime>::Random(rows, rows);
  26. VectorType v1 = VectorType::Random(rows),
  27. vzero = VectorType::Zero(rows);
  28. SquareMatrixType sm1 = SquareMatrixType::Random(rows,rows), sm2(rows,rows);
  29. Scalar x = 0;
  30. while(x == Scalar(0)) x = internal::random<Scalar>();
  31. Index r = internal::random<Index>(0, rows-1),
  32. c = internal::random<Index>(0, cols-1);
  33. m1.coeffRef(r,c) = x;
  34. VERIFY_IS_APPROX(x, m1.coeff(r,c));
  35. m1(r,c) = x;
  36. VERIFY_IS_APPROX(x, m1(r,c));
  37. v1.coeffRef(r) = x;
  38. VERIFY_IS_APPROX(x, v1.coeff(r));
  39. v1(r) = x;
  40. VERIFY_IS_APPROX(x, v1(r));
  41. v1[r] = x;
  42. VERIFY_IS_APPROX(x, v1[r]);
  43. // test fetching with various index types.
  44. Index r1 = internal::random<Index>(0, numext::mini(Index(127),rows-1));
  45. x = v1(static_cast<char>(r1));
  46. x = v1(static_cast<signed char>(r1));
  47. x = v1(static_cast<unsigned char>(r1));
  48. x = v1(static_cast<signed short>(r1));
  49. x = v1(static_cast<unsigned short>(r1));
  50. x = v1(static_cast<signed int>(r1));
  51. x = v1(static_cast<unsigned int>(r1));
  52. x = v1(static_cast<signed long>(r1));
  53. x = v1(static_cast<unsigned long>(r1));
  54. #if EIGEN_HAS_CXX11
  55. x = v1(static_cast<long long int>(r1));
  56. x = v1(static_cast<unsigned long long int>(r1));
  57. #endif
  58. VERIFY_IS_APPROX( v1, v1);
  59. VERIFY_IS_NOT_APPROX( v1, 2*v1);
  60. VERIFY_IS_MUCH_SMALLER_THAN( vzero, v1);
  61. VERIFY_IS_MUCH_SMALLER_THAN( vzero, v1.squaredNorm());
  62. VERIFY_IS_NOT_MUCH_SMALLER_THAN(v1, v1);
  63. VERIFY_IS_APPROX( vzero, v1-v1);
  64. VERIFY_IS_APPROX( m1, m1);
  65. VERIFY_IS_NOT_APPROX( m1, 2*m1);
  66. VERIFY_IS_MUCH_SMALLER_THAN( mzero, m1);
  67. VERIFY_IS_NOT_MUCH_SMALLER_THAN(m1, m1);
  68. VERIFY_IS_APPROX( mzero, m1-m1);
  69. // always test operator() on each read-only expression class,
  70. // in order to check const-qualifiers.
  71. // indeed, if an expression class (here Zero) is meant to be read-only,
  72. // hence has no _write() method, the corresponding MatrixBase method (here zero())
  73. // should return a const-qualified object so that it is the const-qualified
  74. // operator() that gets called, which in turn calls _read().
  75. VERIFY_IS_MUCH_SMALLER_THAN(MatrixType::Zero(rows,cols)(r,c), static_cast<Scalar>(1));
  76. // now test copying a row-vector into a (column-)vector and conversely.
  77. square.col(r) = square.row(r).eval();
  78. Matrix<Scalar, 1, MatrixType::RowsAtCompileTime> rv(rows);
  79. Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> cv(rows);
  80. rv = square.row(r);
  81. cv = square.col(r);
  82. VERIFY_IS_APPROX(rv, cv.transpose());
  83. if(cols!=1 && rows!=1 && MatrixType::SizeAtCompileTime!=Dynamic)
  84. {
  85. VERIFY_RAISES_ASSERT(m1 = (m2.block(0,0, rows-1, cols-1)));
  86. }
  87. if(cols!=1 && rows!=1)
  88. {
  89. VERIFY_RAISES_ASSERT(m1[0]);
  90. VERIFY_RAISES_ASSERT((m1+m1)[0]);
  91. }
  92. VERIFY_IS_APPROX(m3 = m1,m1);
  93. MatrixType m4;
  94. VERIFY_IS_APPROX(m4 = m1,m1);
  95. m3.real() = m1.real();
  96. VERIFY_IS_APPROX(static_cast<const MatrixType&>(m3).real(), static_cast<const MatrixType&>(m1).real());
  97. VERIFY_IS_APPROX(static_cast<const MatrixType&>(m3).real(), m1.real());
  98. // check == / != operators
  99. VERIFY(m1==m1);
  100. VERIFY(m1!=m2);
  101. VERIFY(!(m1==m2));
  102. VERIFY(!(m1!=m1));
  103. m1 = m2;
  104. VERIFY(m1==m2);
  105. VERIFY(!(m1!=m2));
  106. // check automatic transposition
  107. sm2.setZero();
  108. for(Index i=0;i<rows;++i)
  109. sm2.col(i) = sm1.row(i);
  110. VERIFY_IS_APPROX(sm2,sm1.transpose());
  111. sm2.setZero();
  112. for(Index i=0;i<rows;++i)
  113. sm2.col(i).noalias() = sm1.row(i);
  114. VERIFY_IS_APPROX(sm2,sm1.transpose());
  115. sm2.setZero();
  116. for(Index i=0;i<rows;++i)
  117. sm2.col(i).noalias() += sm1.row(i);
  118. VERIFY_IS_APPROX(sm2,sm1.transpose());
  119. sm2.setZero();
  120. for(Index i=0;i<rows;++i)
  121. sm2.col(i).noalias() -= sm1.row(i);
  122. VERIFY_IS_APPROX(sm2,-sm1.transpose());
  123. // check ternary usage
  124. {
  125. bool b = internal::random<int>(0,10)>5;
  126. m3 = b ? m1 : m2;
  127. if(b) VERIFY_IS_APPROX(m3,m1);
  128. else VERIFY_IS_APPROX(m3,m2);
  129. m3 = b ? -m1 : m2;
  130. if(b) VERIFY_IS_APPROX(m3,-m1);
  131. else VERIFY_IS_APPROX(m3,m2);
  132. m3 = b ? m1 : -m2;
  133. if(b) VERIFY_IS_APPROX(m3,m1);
  134. else VERIFY_IS_APPROX(m3,-m2);
  135. }
  136. }
  137. template<typename MatrixType> void basicStuffComplex(const MatrixType& m)
  138. {
  139. typedef typename MatrixType::Scalar Scalar;
  140. typedef typename NumTraits<Scalar>::Real RealScalar;
  141. typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime> RealMatrixType;
  142. Index rows = m.rows();
  143. Index cols = m.cols();
  144. Scalar s1 = internal::random<Scalar>(),
  145. s2 = internal::random<Scalar>();
  146. VERIFY(numext::real(s1)==numext::real_ref(s1));
  147. VERIFY(numext::imag(s1)==numext::imag_ref(s1));
  148. numext::real_ref(s1) = numext::real(s2);
  149. numext::imag_ref(s1) = numext::imag(s2);
  150. VERIFY(internal::isApprox(s1, s2, NumTraits<RealScalar>::epsilon()));
  151. // extended precision in Intel FPUs means that s1 == s2 in the line above is not guaranteed.
  152. RealMatrixType rm1 = RealMatrixType::Random(rows,cols),
  153. rm2 = RealMatrixType::Random(rows,cols);
  154. MatrixType cm(rows,cols);
  155. cm.real() = rm1;
  156. cm.imag() = rm2;
  157. VERIFY_IS_APPROX(static_cast<const MatrixType&>(cm).real(), rm1);
  158. VERIFY_IS_APPROX(static_cast<const MatrixType&>(cm).imag(), rm2);
  159. rm1.setZero();
  160. rm2.setZero();
  161. rm1 = cm.real();
  162. rm2 = cm.imag();
  163. VERIFY_IS_APPROX(static_cast<const MatrixType&>(cm).real(), rm1);
  164. VERIFY_IS_APPROX(static_cast<const MatrixType&>(cm).imag(), rm2);
  165. cm.real().setZero();
  166. VERIFY(static_cast<const MatrixType&>(cm).real().isZero());
  167. VERIFY(!static_cast<const MatrixType&>(cm).imag().isZero());
  168. }
  169. template<typename SrcScalar, typename TgtScalar>
  170. struct casting_test {
  171. static void run() {
  172. Matrix<SrcScalar,4,4> m;
  173. for (int i=0; i<m.rows(); ++i) {
  174. for (int j=0; j<m.cols(); ++j) {
  175. m(i, j) = internal::random_without_cast_overflow<SrcScalar,TgtScalar>::value();
  176. }
  177. }
  178. Matrix<TgtScalar,4,4> n = m.template cast<TgtScalar>();
  179. for (int i=0; i<m.rows(); ++i) {
  180. for (int j=0; j<m.cols(); ++j) {
  181. VERIFY_IS_APPROX(n(i, j), (internal::cast<SrcScalar,TgtScalar>(m(i, j))));
  182. }
  183. }
  184. }
  185. };
  186. template<typename SrcScalar, typename EnableIf = void>
  187. struct casting_test_runner {
  188. static void run() {
  189. casting_test<SrcScalar, bool>::run();
  190. casting_test<SrcScalar, int8_t>::run();
  191. casting_test<SrcScalar, uint8_t>::run();
  192. casting_test<SrcScalar, int16_t>::run();
  193. casting_test<SrcScalar, uint16_t>::run();
  194. casting_test<SrcScalar, int32_t>::run();
  195. casting_test<SrcScalar, uint32_t>::run();
  196. #if EIGEN_HAS_CXX11
  197. casting_test<SrcScalar, int64_t>::run();
  198. casting_test<SrcScalar, uint64_t>::run();
  199. #endif
  200. casting_test<SrcScalar, half>::run();
  201. casting_test<SrcScalar, bfloat16>::run();
  202. casting_test<SrcScalar, float>::run();
  203. casting_test<SrcScalar, double>::run();
  204. casting_test<SrcScalar, std::complex<float> >::run();
  205. casting_test<SrcScalar, std::complex<double> >::run();
  206. }
  207. };
  208. template<typename SrcScalar>
  209. struct casting_test_runner<SrcScalar, typename internal::enable_if<(NumTraits<SrcScalar>::IsComplex)>::type>
  210. {
  211. static void run() {
  212. // Only a few casts from std::complex<T> are defined.
  213. casting_test<SrcScalar, half>::run();
  214. casting_test<SrcScalar, bfloat16>::run();
  215. casting_test<SrcScalar, std::complex<float> >::run();
  216. casting_test<SrcScalar, std::complex<double> >::run();
  217. }
  218. };
  219. void casting_all() {
  220. casting_test_runner<bool>::run();
  221. casting_test_runner<int8_t>::run();
  222. casting_test_runner<uint8_t>::run();
  223. casting_test_runner<int16_t>::run();
  224. casting_test_runner<uint16_t>::run();
  225. casting_test_runner<int32_t>::run();
  226. casting_test_runner<uint32_t>::run();
  227. #if EIGEN_HAS_CXX11
  228. casting_test_runner<int64_t>::run();
  229. casting_test_runner<uint64_t>::run();
  230. #endif
  231. casting_test_runner<half>::run();
  232. casting_test_runner<bfloat16>::run();
  233. casting_test_runner<float>::run();
  234. casting_test_runner<double>::run();
  235. casting_test_runner<std::complex<float> >::run();
  236. casting_test_runner<std::complex<double> >::run();
  237. }
  238. template <typename Scalar>
  239. void fixedSizeMatrixConstruction()
  240. {
  241. Scalar raw[4];
  242. for(int k=0; k<4; ++k)
  243. raw[k] = internal::random<Scalar>();
  244. {
  245. Matrix<Scalar,4,1> m(raw);
  246. Array<Scalar,4,1> a(raw);
  247. for(int k=0; k<4; ++k) VERIFY(m(k) == raw[k]);
  248. for(int k=0; k<4; ++k) VERIFY(a(k) == raw[k]);
  249. VERIFY_IS_EQUAL(m,(Matrix<Scalar,4,1>(raw[0],raw[1],raw[2],raw[3])));
  250. VERIFY((a==(Array<Scalar,4,1>(raw[0],raw[1],raw[2],raw[3]))).all());
  251. }
  252. {
  253. Matrix<Scalar,3,1> m(raw);
  254. Array<Scalar,3,1> a(raw);
  255. for(int k=0; k<3; ++k) VERIFY(m(k) == raw[k]);
  256. for(int k=0; k<3; ++k) VERIFY(a(k) == raw[k]);
  257. VERIFY_IS_EQUAL(m,(Matrix<Scalar,3,1>(raw[0],raw[1],raw[2])));
  258. VERIFY((a==Array<Scalar,3,1>(raw[0],raw[1],raw[2])).all());
  259. }
  260. {
  261. Matrix<Scalar,2,1> m(raw), m2( (DenseIndex(raw[0])), (DenseIndex(raw[1])) );
  262. Array<Scalar,2,1> a(raw), a2( (DenseIndex(raw[0])), (DenseIndex(raw[1])) );
  263. for(int k=0; k<2; ++k) VERIFY(m(k) == raw[k]);
  264. for(int k=0; k<2; ++k) VERIFY(a(k) == raw[k]);
  265. VERIFY_IS_EQUAL(m,(Matrix<Scalar,2,1>(raw[0],raw[1])));
  266. VERIFY((a==Array<Scalar,2,1>(raw[0],raw[1])).all());
  267. for(int k=0; k<2; ++k) VERIFY(m2(k) == DenseIndex(raw[k]));
  268. for(int k=0; k<2; ++k) VERIFY(a2(k) == DenseIndex(raw[k]));
  269. }
  270. {
  271. Matrix<Scalar,1,2> m(raw),
  272. m2( (DenseIndex(raw[0])), (DenseIndex(raw[1])) ),
  273. m3( (int(raw[0])), (int(raw[1])) ),
  274. m4( (float(raw[0])), (float(raw[1])) );
  275. Array<Scalar,1,2> a(raw), a2( (DenseIndex(raw[0])), (DenseIndex(raw[1])) );
  276. for(int k=0; k<2; ++k) VERIFY(m(k) == raw[k]);
  277. for(int k=0; k<2; ++k) VERIFY(a(k) == raw[k]);
  278. VERIFY_IS_EQUAL(m,(Matrix<Scalar,1,2>(raw[0],raw[1])));
  279. VERIFY((a==Array<Scalar,1,2>(raw[0],raw[1])).all());
  280. for(int k=0; k<2; ++k) VERIFY(m2(k) == DenseIndex(raw[k]));
  281. for(int k=0; k<2; ++k) VERIFY(a2(k) == DenseIndex(raw[k]));
  282. for(int k=0; k<2; ++k) VERIFY(m3(k) == int(raw[k]));
  283. for(int k=0; k<2; ++k) VERIFY((m4(k)) == Scalar(float(raw[k])));
  284. }
  285. {
  286. Matrix<Scalar,1,1> m(raw), m1(raw[0]), m2( (DenseIndex(raw[0])) ), m3( (int(raw[0])) );
  287. Array<Scalar,1,1> a(raw), a1(raw[0]), a2( (DenseIndex(raw[0])) );
  288. VERIFY(m(0) == raw[0]);
  289. VERIFY(a(0) == raw[0]);
  290. VERIFY(m1(0) == raw[0]);
  291. VERIFY(a1(0) == raw[0]);
  292. VERIFY(m2(0) == DenseIndex(raw[0]));
  293. VERIFY(a2(0) == DenseIndex(raw[0]));
  294. VERIFY(m3(0) == int(raw[0]));
  295. VERIFY_IS_EQUAL(m,(Matrix<Scalar,1,1>(raw[0])));
  296. VERIFY((a==Array<Scalar,1,1>(raw[0])).all());
  297. }
  298. }
  299. EIGEN_DECLARE_TEST(basicstuff)
  300. {
  301. for(int i = 0; i < g_repeat; i++) {
  302. CALL_SUBTEST_1( basicStuff(Matrix<float, 1, 1>()) );
  303. CALL_SUBTEST_2( basicStuff(Matrix4d()) );
  304. CALL_SUBTEST_3( basicStuff(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  305. CALL_SUBTEST_4( basicStuff(MatrixXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  306. CALL_SUBTEST_5( basicStuff(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  307. CALL_SUBTEST_6( basicStuff(Matrix<float, 100, 100>()) );
  308. CALL_SUBTEST_7( basicStuff(Matrix<long double,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE),internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  309. CALL_SUBTEST_8( casting_all() );
  310. CALL_SUBTEST_3( basicStuffComplex(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  311. CALL_SUBTEST_5( basicStuffComplex(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  312. }
  313. CALL_SUBTEST_1(fixedSizeMatrixConstruction<unsigned char>());
  314. CALL_SUBTEST_1(fixedSizeMatrixConstruction<float>());
  315. CALL_SUBTEST_1(fixedSizeMatrixConstruction<double>());
  316. CALL_SUBTEST_1(fixedSizeMatrixConstruction<int>());
  317. CALL_SUBTEST_1(fixedSizeMatrixConstruction<long int>());
  318. CALL_SUBTEST_1(fixedSizeMatrixConstruction<std::ptrdiff_t>());
  319. }