initializer_list_construction.cpp 12 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385
  1. // This file is part of Eigen, a lightweight C++ template library
  2. // for linear algebra.
  3. //
  4. // Copyright (C) 2019 David Tellenbach <david.tellenbach@tellnotes.org>
  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. template<typename Scalar, bool is_integer = NumTraits<Scalar>::IsInteger>
  12. struct TestMethodDispatching {
  13. static void run() {}
  14. };
  15. template<typename Scalar>
  16. struct TestMethodDispatching<Scalar, 1> {
  17. static void run()
  18. {
  19. {
  20. Matrix<Scalar, Dynamic, Dynamic> m {3, 4};
  21. Array<Scalar, Dynamic, Dynamic> a {3, 4};
  22. VERIFY(m.rows() == 3);
  23. VERIFY(m.cols() == 4);
  24. VERIFY(a.rows() == 3);
  25. VERIFY(a.cols() == 4);
  26. }
  27. {
  28. Matrix<Scalar, 1, 2> m {3, 4};
  29. Array<Scalar, 1, 2> a {3, 4};
  30. VERIFY(m(0) == 3);
  31. VERIFY(m(1) == 4);
  32. VERIFY(a(0) == 3);
  33. VERIFY(a(1) == 4);
  34. }
  35. {
  36. Matrix<Scalar, 2, 1> m {3, 4};
  37. Array<Scalar, 2, 1> a {3, 4};
  38. VERIFY(m(0) == 3);
  39. VERIFY(m(1) == 4);
  40. VERIFY(a(0) == 3);
  41. VERIFY(a(1) == 4);
  42. }
  43. }
  44. };
  45. template<typename Vec4, typename Vec5> void fixedsizeVariadicVectorConstruction2()
  46. {
  47. {
  48. Vec4 ref = Vec4::Random();
  49. Vec4 v{ ref[0], ref[1], ref[2], ref[3] };
  50. VERIFY_IS_APPROX(v, ref);
  51. VERIFY_IS_APPROX(v, (Vec4( ref[0], ref[1], ref[2], ref[3] )));
  52. VERIFY_IS_APPROX(v, (Vec4({ref[0], ref[1], ref[2], ref[3]})));
  53. Vec4 v2 = { ref[0], ref[1], ref[2], ref[3] };
  54. VERIFY_IS_APPROX(v2, ref);
  55. }
  56. {
  57. Vec5 ref = Vec5::Random();
  58. Vec5 v{ ref[0], ref[1], ref[2], ref[3], ref[4] };
  59. VERIFY_IS_APPROX(v, ref);
  60. VERIFY_IS_APPROX(v, (Vec5( ref[0], ref[1], ref[2], ref[3], ref[4] )));
  61. VERIFY_IS_APPROX(v, (Vec5({ref[0], ref[1], ref[2], ref[3], ref[4]})));
  62. Vec5 v2 = { ref[0], ref[1], ref[2], ref[3], ref[4] };
  63. VERIFY_IS_APPROX(v2, ref);
  64. }
  65. }
  66. #define CHECK_MIXSCALAR_V5_APPROX(V, A0, A1, A2, A3, A4) { \
  67. VERIFY_IS_APPROX(V[0], Scalar(A0) ); \
  68. VERIFY_IS_APPROX(V[1], Scalar(A1) ); \
  69. VERIFY_IS_APPROX(V[2], Scalar(A2) ); \
  70. VERIFY_IS_APPROX(V[3], Scalar(A3) ); \
  71. VERIFY_IS_APPROX(V[4], Scalar(A4) ); \
  72. }
  73. #define CHECK_MIXSCALAR_V5(VEC5, A0, A1, A2, A3, A4) { \
  74. typedef VEC5::Scalar Scalar; \
  75. VEC5 v = { A0 , A1 , A2 , A3 , A4 }; \
  76. CHECK_MIXSCALAR_V5_APPROX(v, A0 , A1 , A2 , A3 , A4); \
  77. }
  78. template<int> void fixedsizeVariadicVectorConstruction3()
  79. {
  80. typedef Matrix<double,5,1> Vec5;
  81. typedef Array<float,5,1> Arr5;
  82. CHECK_MIXSCALAR_V5(Vec5, 1, 2., -3, 4.121, 5.53252);
  83. CHECK_MIXSCALAR_V5(Arr5, 1, 2., 3.12f, 4.121, 5.53252);
  84. }
  85. template<typename Scalar> void fixedsizeVariadicVectorConstruction()
  86. {
  87. CALL_SUBTEST(( fixedsizeVariadicVectorConstruction2<Matrix<Scalar,4,1>, Matrix<Scalar,5,1> >() ));
  88. CALL_SUBTEST(( fixedsizeVariadicVectorConstruction2<Matrix<Scalar,1,4>, Matrix<Scalar,1,5> >() ));
  89. CALL_SUBTEST(( fixedsizeVariadicVectorConstruction2<Array<Scalar,4,1>, Array<Scalar,5,1> >() ));
  90. CALL_SUBTEST(( fixedsizeVariadicVectorConstruction2<Array<Scalar,1,4>, Array<Scalar,1,5> >() ));
  91. }
  92. template<typename Scalar> void initializerListVectorConstruction()
  93. {
  94. Scalar raw[4];
  95. for(int k = 0; k < 4; ++k) {
  96. raw[k] = internal::random<Scalar>();
  97. }
  98. {
  99. Matrix<Scalar, 4, 1> m { {raw[0]}, {raw[1]},{raw[2]},{raw[3]} };
  100. Array<Scalar, 4, 1> a { {raw[0]}, {raw[1]}, {raw[2]}, {raw[3]} };
  101. for(int k = 0; k < 4; ++k) {
  102. VERIFY(m(k) == raw[k]);
  103. }
  104. for(int k = 0; k < 4; ++k) {
  105. VERIFY(a(k) == raw[k]);
  106. }
  107. VERIFY_IS_EQUAL(m, (Matrix<Scalar,4,1>({ {raw[0]}, {raw[1]}, {raw[2]}, {raw[3]} })));
  108. VERIFY((a == (Array<Scalar,4,1>({ {raw[0]}, {raw[1]}, {raw[2]}, {raw[3]} }))).all());
  109. }
  110. {
  111. Matrix<Scalar, 1, 4> m { {raw[0], raw[1], raw[2], raw[3]} };
  112. Array<Scalar, 1, 4> a { {raw[0], raw[1], raw[2], raw[3]} };
  113. for(int k = 0; k < 4; ++k) {
  114. VERIFY(m(k) == raw[k]);
  115. }
  116. for(int k = 0; k < 4; ++k) {
  117. VERIFY(a(k) == raw[k]);
  118. }
  119. VERIFY_IS_EQUAL(m, (Matrix<Scalar, 1, 4>({{raw[0],raw[1],raw[2],raw[3]}})));
  120. VERIFY((a == (Array<Scalar, 1, 4>({{raw[0],raw[1],raw[2],raw[3]}}))).all());
  121. }
  122. {
  123. Matrix<Scalar, 4, Dynamic> m { {raw[0]}, {raw[1]}, {raw[2]}, {raw[3]} };
  124. Array<Scalar, 4, Dynamic> a { {raw[0]}, {raw[1]}, {raw[2]}, {raw[3]} };
  125. for(int k=0; k < 4; ++k) {
  126. VERIFY(m(k) == raw[k]);
  127. }
  128. for(int k=0; k < 4; ++k) {
  129. VERIFY(a(k) == raw[k]);
  130. }
  131. VERIFY_IS_EQUAL(m, (Matrix<Scalar, 4, Dynamic>({ {raw[0]}, {raw[1]}, {raw[2]}, {raw[3]} })));
  132. VERIFY((a == (Array<Scalar, 4, Dynamic>({ {raw[0]}, {raw[1]}, {raw[2]}, {raw[3]} }))).all());
  133. }
  134. {
  135. Matrix<Scalar, Dynamic, 4> m {{raw[0],raw[1],raw[2],raw[3]}};
  136. Array<Scalar, Dynamic, 4> a {{raw[0],raw[1],raw[2],raw[3]}};
  137. for(int k=0; k < 4; ++k) {
  138. VERIFY(m(k) == raw[k]);
  139. }
  140. for(int k=0; k < 4; ++k) {
  141. VERIFY(a(k) == raw[k]);
  142. }
  143. VERIFY_IS_EQUAL(m, (Matrix<Scalar, Dynamic, 4>({{raw[0],raw[1],raw[2],raw[3]}})));
  144. VERIFY((a == (Array<Scalar, Dynamic, 4>({{raw[0],raw[1],raw[2],raw[3]}}))).all());
  145. }
  146. }
  147. template<typename Scalar> void initializerListMatrixConstruction()
  148. {
  149. const Index RowsAtCompileTime = 5;
  150. const Index ColsAtCompileTime = 4;
  151. const Index SizeAtCompileTime = RowsAtCompileTime * ColsAtCompileTime;
  152. Scalar raw[SizeAtCompileTime];
  153. for (int i = 0; i < SizeAtCompileTime; ++i) {
  154. raw[i] = internal::random<Scalar>();
  155. }
  156. {
  157. Matrix<Scalar, Dynamic, Dynamic> m {};
  158. VERIFY(m.cols() == 0);
  159. VERIFY(m.rows() == 0);
  160. VERIFY_IS_EQUAL(m, (Matrix<Scalar, Dynamic, Dynamic>()));
  161. }
  162. {
  163. Matrix<Scalar, 5, 4> m {
  164. {raw[0], raw[1], raw[2], raw[3]},
  165. {raw[4], raw[5], raw[6], raw[7]},
  166. {raw[8], raw[9], raw[10], raw[11]},
  167. {raw[12], raw[13], raw[14], raw[15]},
  168. {raw[16], raw[17], raw[18], raw[19]}
  169. };
  170. Matrix<Scalar, 5, 4> m2;
  171. m2 << raw[0], raw[1], raw[2], raw[3],
  172. raw[4], raw[5], raw[6], raw[7],
  173. raw[8], raw[9], raw[10], raw[11],
  174. raw[12], raw[13], raw[14], raw[15],
  175. raw[16], raw[17], raw[18], raw[19];
  176. int k = 0;
  177. for(int i = 0; i < RowsAtCompileTime; ++i) {
  178. for (int j = 0; j < ColsAtCompileTime; ++j) {
  179. VERIFY(m(i, j) == raw[k]);
  180. ++k;
  181. }
  182. }
  183. VERIFY_IS_EQUAL(m, m2);
  184. }
  185. {
  186. Matrix<Scalar, Dynamic, Dynamic> m{
  187. {raw[0], raw[1], raw[2], raw[3]},
  188. {raw[4], raw[5], raw[6], raw[7]},
  189. {raw[8], raw[9], raw[10], raw[11]},
  190. {raw[12], raw[13], raw[14], raw[15]},
  191. {raw[16], raw[17], raw[18], raw[19]}
  192. };
  193. VERIFY(m.cols() == 4);
  194. VERIFY(m.rows() == 5);
  195. int k = 0;
  196. for(int i = 0; i < RowsAtCompileTime; ++i) {
  197. for (int j = 0; j < ColsAtCompileTime; ++j) {
  198. VERIFY(m(i, j) == raw[k]);
  199. ++k;
  200. }
  201. }
  202. Matrix<Scalar, Dynamic, Dynamic> m2(RowsAtCompileTime, ColsAtCompileTime);
  203. k = 0;
  204. for(int i = 0; i < RowsAtCompileTime; ++i) {
  205. for (int j = 0; j < ColsAtCompileTime; ++j) {
  206. m2(i, j) = raw[k];
  207. ++k;
  208. }
  209. }
  210. VERIFY_IS_EQUAL(m, m2);
  211. }
  212. }
  213. template<typename Scalar> void initializerListArrayConstruction()
  214. {
  215. const Index RowsAtCompileTime = 5;
  216. const Index ColsAtCompileTime = 4;
  217. const Index SizeAtCompileTime = RowsAtCompileTime * ColsAtCompileTime;
  218. Scalar raw[SizeAtCompileTime];
  219. for (int i = 0; i < SizeAtCompileTime; ++i) {
  220. raw[i] = internal::random<Scalar>();
  221. }
  222. {
  223. Array<Scalar, Dynamic, Dynamic> a {};
  224. VERIFY(a.cols() == 0);
  225. VERIFY(a.rows() == 0);
  226. }
  227. {
  228. Array<Scalar, 5, 4> m {
  229. {raw[0], raw[1], raw[2], raw[3]},
  230. {raw[4], raw[5], raw[6], raw[7]},
  231. {raw[8], raw[9], raw[10], raw[11]},
  232. {raw[12], raw[13], raw[14], raw[15]},
  233. {raw[16], raw[17], raw[18], raw[19]}
  234. };
  235. Array<Scalar, 5, 4> m2;
  236. m2 << raw[0], raw[1], raw[2], raw[3],
  237. raw[4], raw[5], raw[6], raw[7],
  238. raw[8], raw[9], raw[10], raw[11],
  239. raw[12], raw[13], raw[14], raw[15],
  240. raw[16], raw[17], raw[18], raw[19];
  241. int k = 0;
  242. for(int i = 0; i < RowsAtCompileTime; ++i) {
  243. for (int j = 0; j < ColsAtCompileTime; ++j) {
  244. VERIFY(m(i, j) == raw[k]);
  245. ++k;
  246. }
  247. }
  248. VERIFY_IS_APPROX(m, m2);
  249. }
  250. {
  251. Array<Scalar, Dynamic, Dynamic> m {
  252. {raw[0], raw[1], raw[2], raw[3]},
  253. {raw[4], raw[5], raw[6], raw[7]},
  254. {raw[8], raw[9], raw[10], raw[11]},
  255. {raw[12], raw[13], raw[14], raw[15]},
  256. {raw[16], raw[17], raw[18], raw[19]}
  257. };
  258. VERIFY(m.cols() == 4);
  259. VERIFY(m.rows() == 5);
  260. int k = 0;
  261. for(int i = 0; i < RowsAtCompileTime; ++i) {
  262. for (int j = 0; j < ColsAtCompileTime; ++j) {
  263. VERIFY(m(i, j) == raw[k]);
  264. ++k;
  265. }
  266. }
  267. Array<Scalar, Dynamic, Dynamic> m2(RowsAtCompileTime, ColsAtCompileTime);
  268. k = 0;
  269. for(int i = 0; i < RowsAtCompileTime; ++i) {
  270. for (int j = 0; j < ColsAtCompileTime; ++j) {
  271. m2(i, j) = raw[k];
  272. ++k;
  273. }
  274. }
  275. VERIFY_IS_APPROX(m, m2);
  276. }
  277. }
  278. template<typename Scalar> void dynamicVectorConstruction()
  279. {
  280. const Index size = 4;
  281. Scalar raw[size];
  282. for (int i = 0; i < size; ++i) {
  283. raw[i] = internal::random<Scalar>();
  284. }
  285. typedef Matrix<Scalar, Dynamic, 1> VectorX;
  286. {
  287. VectorX v {{raw[0], raw[1], raw[2], raw[3]}};
  288. for (int i = 0; i < size; ++i) {
  289. VERIFY(v(i) == raw[i]);
  290. }
  291. VERIFY(v.rows() == size);
  292. VERIFY(v.cols() == 1);
  293. VERIFY_IS_EQUAL(v, (VectorX {{raw[0], raw[1], raw[2], raw[3]}}));
  294. }
  295. {
  296. VERIFY_RAISES_ASSERT((VectorX {raw[0], raw[1], raw[2], raw[3]}));
  297. }
  298. {
  299. VERIFY_RAISES_ASSERT((VectorX {
  300. {raw[0], raw[1], raw[2], raw[3]},
  301. {raw[0], raw[1], raw[2], raw[3]},
  302. }));
  303. }
  304. }
  305. EIGEN_DECLARE_TEST(initializer_list_construction)
  306. {
  307. CALL_SUBTEST_1(initializerListVectorConstruction<unsigned char>());
  308. CALL_SUBTEST_1(initializerListVectorConstruction<float>());
  309. CALL_SUBTEST_1(initializerListVectorConstruction<double>());
  310. CALL_SUBTEST_1(initializerListVectorConstruction<int>());
  311. CALL_SUBTEST_1(initializerListVectorConstruction<long int>());
  312. CALL_SUBTEST_1(initializerListVectorConstruction<std::ptrdiff_t>());
  313. CALL_SUBTEST_1(initializerListVectorConstruction<std::complex<double>>());
  314. CALL_SUBTEST_1(initializerListVectorConstruction<std::complex<float>>());
  315. CALL_SUBTEST_2(initializerListMatrixConstruction<unsigned char>());
  316. CALL_SUBTEST_2(initializerListMatrixConstruction<float>());
  317. CALL_SUBTEST_2(initializerListMatrixConstruction<double>());
  318. CALL_SUBTEST_2(initializerListMatrixConstruction<int>());
  319. CALL_SUBTEST_2(initializerListMatrixConstruction<long int>());
  320. CALL_SUBTEST_2(initializerListMatrixConstruction<std::ptrdiff_t>());
  321. CALL_SUBTEST_2(initializerListMatrixConstruction<std::complex<double>>());
  322. CALL_SUBTEST_2(initializerListMatrixConstruction<std::complex<float>>());
  323. CALL_SUBTEST_3(initializerListArrayConstruction<unsigned char>());
  324. CALL_SUBTEST_3(initializerListArrayConstruction<float>());
  325. CALL_SUBTEST_3(initializerListArrayConstruction<double>());
  326. CALL_SUBTEST_3(initializerListArrayConstruction<int>());
  327. CALL_SUBTEST_3(initializerListArrayConstruction<long int>());
  328. CALL_SUBTEST_3(initializerListArrayConstruction<std::ptrdiff_t>());
  329. CALL_SUBTEST_3(initializerListArrayConstruction<std::complex<double>>());
  330. CALL_SUBTEST_3(initializerListArrayConstruction<std::complex<float>>());
  331. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction<unsigned char>());
  332. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction<float>());
  333. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction<double>());
  334. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction<int>());
  335. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction<long int>());
  336. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction<std::ptrdiff_t>());
  337. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction<std::complex<double>>());
  338. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction<std::complex<float>>());
  339. CALL_SUBTEST_4(fixedsizeVariadicVectorConstruction3<0>());
  340. CALL_SUBTEST_5(TestMethodDispatching<int>::run());
  341. CALL_SUBTEST_5(TestMethodDispatching<long int>::run());
  342. CALL_SUBTEST_6(dynamicVectorConstruction<unsigned char>());
  343. CALL_SUBTEST_6(dynamicVectorConstruction<float>());
  344. CALL_SUBTEST_6(dynamicVectorConstruction<double>());
  345. CALL_SUBTEST_6(dynamicVectorConstruction<int>());
  346. CALL_SUBTEST_6(dynamicVectorConstruction<long int>());
  347. CALL_SUBTEST_6(dynamicVectorConstruction<std::ptrdiff_t>());
  348. CALL_SUBTEST_6(dynamicVectorConstruction<std::complex<double>>());
  349. CALL_SUBTEST_6(dynamicVectorConstruction<std::complex<float>>());
  350. }