array_cwise.cpp 27 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576577578579580581582583584585586587588589590591592593594595596597598599600601602603604605606607608609610611612613614615616617618619620621622623624625626627628629630631632633634635636637638639640641642643644645646647648649650651652653654655656657658659660661662663664665666667668669670671672673674675676677678679680681682683684685686687688689690691692693694695696697698699700701702703704705706707708709710
  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. // Test the corner cases of pow(x, y) for real types.
  11. template<typename Scalar>
  12. void pow_test() {
  13. const Scalar zero = Scalar(0);
  14. const Scalar eps = Eigen::NumTraits<Scalar>::epsilon();
  15. const Scalar one = Scalar(1);
  16. const Scalar two = Scalar(2);
  17. const Scalar three = Scalar(3);
  18. const Scalar sqrt_half = Scalar(std::sqrt(0.5));
  19. const Scalar sqrt2 = Scalar(std::sqrt(2));
  20. const Scalar inf = Eigen::NumTraits<Scalar>::infinity();
  21. const Scalar nan = Eigen::NumTraits<Scalar>::quiet_NaN();
  22. const Scalar denorm_min = std::numeric_limits<Scalar>::denorm_min();
  23. const Scalar min = (std::numeric_limits<Scalar>::min)();
  24. const Scalar max = (std::numeric_limits<Scalar>::max)();
  25. const Scalar max_exp = (static_cast<Scalar>(int(Eigen::NumTraits<Scalar>::max_exponent())) * Scalar(EIGEN_LN2)) / eps;
  26. const static Scalar abs_vals[] = {zero,
  27. denorm_min,
  28. min,
  29. eps,
  30. sqrt_half,
  31. one,
  32. sqrt2,
  33. two,
  34. three,
  35. max_exp,
  36. max,
  37. inf,
  38. nan};
  39. const int abs_cases = 13;
  40. const int num_cases = 2*abs_cases * 2*abs_cases;
  41. // Repeat the same value to make sure we hit the vectorized path.
  42. const int num_repeats = 32;
  43. Array<Scalar, Dynamic, Dynamic> x(num_repeats, num_cases);
  44. Array<Scalar, Dynamic, Dynamic> y(num_repeats, num_cases);
  45. int count = 0;
  46. for (int i = 0; i < abs_cases; ++i) {
  47. const Scalar abs_x = abs_vals[i];
  48. for (int sign_x = 0; sign_x < 2; ++sign_x) {
  49. Scalar x_case = sign_x == 0 ? -abs_x : abs_x;
  50. for (int j = 0; j < abs_cases; ++j) {
  51. const Scalar abs_y = abs_vals[j];
  52. for (int sign_y = 0; sign_y < 2; ++sign_y) {
  53. Scalar y_case = sign_y == 0 ? -abs_y : abs_y;
  54. for (int repeat = 0; repeat < num_repeats; ++repeat) {
  55. x(repeat, count) = x_case;
  56. y(repeat, count) = y_case;
  57. }
  58. ++count;
  59. }
  60. }
  61. }
  62. }
  63. Array<Scalar, Dynamic, Dynamic> actual = x.pow(y);
  64. const Scalar tol = test_precision<Scalar>();
  65. bool all_pass = true;
  66. for (int i = 0; i < 1; ++i) {
  67. for (int j = 0; j < num_cases; ++j) {
  68. Scalar e = static_cast<Scalar>(std::pow(x(i,j), y(i,j)));
  69. Scalar a = actual(i, j);
  70. bool fail = !(a==e) && !internal::isApprox(a, e, tol) && !((numext::isnan)(a) && (numext::isnan)(e));
  71. all_pass &= !fail;
  72. if (fail) {
  73. std::cout << "pow(" << x(i,j) << "," << y(i,j) << ") = " << a << " != " << e << std::endl;
  74. }
  75. }
  76. }
  77. VERIFY(all_pass);
  78. }
  79. template<typename ArrayType> void array(const ArrayType& m)
  80. {
  81. typedef typename ArrayType::Scalar Scalar;
  82. typedef typename ArrayType::RealScalar RealScalar;
  83. typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
  84. typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
  85. Index rows = m.rows();
  86. Index cols = m.cols();
  87. ArrayType m1 = ArrayType::Random(rows, cols),
  88. m2 = ArrayType::Random(rows, cols),
  89. m3(rows, cols);
  90. ArrayType m4 = m1; // copy constructor
  91. VERIFY_IS_APPROX(m1, m4);
  92. ColVectorType cv1 = ColVectorType::Random(rows);
  93. RowVectorType rv1 = RowVectorType::Random(cols);
  94. Scalar s1 = internal::random<Scalar>(),
  95. s2 = internal::random<Scalar>();
  96. // scalar addition
  97. VERIFY_IS_APPROX(m1 + s1, s1 + m1);
  98. VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
  99. VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
  100. VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
  101. VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
  102. VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
  103. m3 = m1;
  104. m3 += s2;
  105. VERIFY_IS_APPROX(m3, m1 + s2);
  106. m3 = m1;
  107. m3 -= s1;
  108. VERIFY_IS_APPROX(m3, m1 - s1);
  109. // scalar operators via Maps
  110. m3 = m1;
  111. ArrayType::Map(m1.data(), m1.rows(), m1.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
  112. VERIFY_IS_APPROX(m1, m3 - m2);
  113. m3 = m1;
  114. ArrayType::Map(m1.data(), m1.rows(), m1.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
  115. VERIFY_IS_APPROX(m1, m3 + m2);
  116. m3 = m1;
  117. ArrayType::Map(m1.data(), m1.rows(), m1.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
  118. VERIFY_IS_APPROX(m1, m3 * m2);
  119. m3 = m1;
  120. m2 = ArrayType::Random(rows,cols);
  121. m2 = (m2==0).select(1,m2);
  122. ArrayType::Map(m1.data(), m1.rows(), m1.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
  123. VERIFY_IS_APPROX(m1, m3 / m2);
  124. // reductions
  125. VERIFY_IS_APPROX(m1.abs().colwise().sum().sum(), m1.abs().sum());
  126. VERIFY_IS_APPROX(m1.abs().rowwise().sum().sum(), m1.abs().sum());
  127. using std::abs;
  128. VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.colwise().sum().sum() - m1.sum()), m1.abs().sum());
  129. VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.rowwise().sum().sum() - m1.sum()), m1.abs().sum());
  130. if (!internal::isMuchSmallerThan(abs(m1.sum() - (m1+m2).sum()), m1.abs().sum(), test_precision<Scalar>()))
  131. VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
  132. VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
  133. // vector-wise ops
  134. m3 = m1;
  135. VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
  136. m3 = m1;
  137. VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
  138. m3 = m1;
  139. VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
  140. m3 = m1;
  141. VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
  142. // Conversion from scalar
  143. VERIFY_IS_APPROX((m3 = s1), ArrayType::Constant(rows,cols,s1));
  144. VERIFY_IS_APPROX((m3 = 1), ArrayType::Constant(rows,cols,1));
  145. VERIFY_IS_APPROX((m3.topLeftCorner(rows,cols) = 1), ArrayType::Constant(rows,cols,1));
  146. typedef Array<Scalar,
  147. ArrayType::RowsAtCompileTime==Dynamic?2:ArrayType::RowsAtCompileTime,
  148. ArrayType::ColsAtCompileTime==Dynamic?2:ArrayType::ColsAtCompileTime,
  149. ArrayType::Options> FixedArrayType;
  150. {
  151. FixedArrayType f1(s1);
  152. VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
  153. FixedArrayType f2(numext::real(s1));
  154. VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
  155. FixedArrayType f3((int)100*numext::real(s1));
  156. VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
  157. f1.setRandom();
  158. FixedArrayType f4(f1.data());
  159. VERIFY_IS_APPROX(f4, f1);
  160. }
  161. #if EIGEN_HAS_CXX11
  162. {
  163. FixedArrayType f1{s1};
  164. VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
  165. FixedArrayType f2{numext::real(s1)};
  166. VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
  167. FixedArrayType f3{(int)100*numext::real(s1)};
  168. VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
  169. f1.setRandom();
  170. FixedArrayType f4{f1.data()};
  171. VERIFY_IS_APPROX(f4, f1);
  172. }
  173. #endif
  174. // pow
  175. VERIFY_IS_APPROX(m1.pow(2), m1.square());
  176. VERIFY_IS_APPROX(pow(m1,2), m1.square());
  177. VERIFY_IS_APPROX(m1.pow(3), m1.cube());
  178. VERIFY_IS_APPROX(pow(m1,3), m1.cube());
  179. VERIFY_IS_APPROX((-m1).pow(3), -m1.cube());
  180. VERIFY_IS_APPROX(pow(2*m1,3), 8*m1.cube());
  181. ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2));
  182. VERIFY_IS_APPROX(Eigen::pow(m1,exponents), m1.square());
  183. VERIFY_IS_APPROX(m1.pow(exponents), m1.square());
  184. VERIFY_IS_APPROX(Eigen::pow(2*m1,exponents), 4*m1.square());
  185. VERIFY_IS_APPROX((2*m1).pow(exponents), 4*m1.square());
  186. VERIFY_IS_APPROX(Eigen::pow(m1,2*exponents), m1.square().square());
  187. VERIFY_IS_APPROX(m1.pow(2*exponents), m1.square().square());
  188. VERIFY_IS_APPROX(Eigen::pow(m1(0,0), exponents), ArrayType::Constant(rows,cols,m1(0,0)*m1(0,0)));
  189. // Check possible conflicts with 1D ctor
  190. typedef Array<Scalar, Dynamic, 1> OneDArrayType;
  191. {
  192. OneDArrayType o1(rows);
  193. VERIFY(o1.size()==rows);
  194. OneDArrayType o2(static_cast<int>(rows));
  195. VERIFY(o2.size()==rows);
  196. }
  197. #if EIGEN_HAS_CXX11
  198. {
  199. OneDArrayType o1{rows};
  200. VERIFY(o1.size()==rows);
  201. OneDArrayType o4{int(rows)};
  202. VERIFY(o4.size()==rows);
  203. }
  204. #endif
  205. // Check possible conflicts with 2D ctor
  206. typedef Array<Scalar, Dynamic, Dynamic> TwoDArrayType;
  207. typedef Array<Scalar, 2, 1> ArrayType2;
  208. {
  209. TwoDArrayType o1(rows,cols);
  210. VERIFY(o1.rows()==rows);
  211. VERIFY(o1.cols()==cols);
  212. TwoDArrayType o2(static_cast<int>(rows),static_cast<int>(cols));
  213. VERIFY(o2.rows()==rows);
  214. VERIFY(o2.cols()==cols);
  215. ArrayType2 o3(rows,cols);
  216. VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
  217. ArrayType2 o4(static_cast<int>(rows),static_cast<int>(cols));
  218. VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
  219. }
  220. #if EIGEN_HAS_CXX11
  221. {
  222. TwoDArrayType o1{rows,cols};
  223. VERIFY(o1.rows()==rows);
  224. VERIFY(o1.cols()==cols);
  225. TwoDArrayType o2{int(rows),int(cols)};
  226. VERIFY(o2.rows()==rows);
  227. VERIFY(o2.cols()==cols);
  228. ArrayType2 o3{rows,cols};
  229. VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
  230. ArrayType2 o4{int(rows),int(cols)};
  231. VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
  232. }
  233. #endif
  234. }
  235. template<typename ArrayType> void comparisons(const ArrayType& m)
  236. {
  237. using std::abs;
  238. typedef typename ArrayType::Scalar Scalar;
  239. typedef typename NumTraits<Scalar>::Real RealScalar;
  240. Index rows = m.rows();
  241. Index cols = m.cols();
  242. Index r = internal::random<Index>(0, rows-1),
  243. c = internal::random<Index>(0, cols-1);
  244. ArrayType m1 = ArrayType::Random(rows, cols),
  245. m2 = ArrayType::Random(rows, cols),
  246. m3(rows, cols),
  247. m4 = m1;
  248. m4 = (m4.abs()==Scalar(0)).select(1,m4);
  249. VERIFY(((m1 + Scalar(1)) > m1).all());
  250. VERIFY(((m1 - Scalar(1)) < m1).all());
  251. if (rows*cols>1)
  252. {
  253. m3 = m1;
  254. m3(r,c) += 1;
  255. VERIFY(! (m1 < m3).all() );
  256. VERIFY(! (m1 > m3).all() );
  257. }
  258. VERIFY(!(m1 > m2 && m1 < m2).any());
  259. VERIFY((m1 <= m2 || m1 >= m2).all());
  260. // comparisons array to scalar
  261. VERIFY( (m1 != (m1(r,c)+1) ).any() );
  262. VERIFY( (m1 > (m1(r,c)-1) ).any() );
  263. VERIFY( (m1 < (m1(r,c)+1) ).any() );
  264. VERIFY( (m1 == m1(r,c) ).any() );
  265. // comparisons scalar to array
  266. VERIFY( ( (m1(r,c)+1) != m1).any() );
  267. VERIFY( ( (m1(r,c)-1) < m1).any() );
  268. VERIFY( ( (m1(r,c)+1) > m1).any() );
  269. VERIFY( ( m1(r,c) == m1).any() );
  270. // test Select
  271. VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
  272. VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
  273. Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
  274. for (int j=0; j<cols; ++j)
  275. for (int i=0; i<rows; ++i)
  276. m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
  277. VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
  278. .select(ArrayType::Zero(rows,cols),m1), m3);
  279. // shorter versions:
  280. VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
  281. .select(0,m1), m3);
  282. VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
  283. .select(m1,0), m3);
  284. // even shorter version:
  285. VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
  286. // count
  287. VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
  288. // and/or
  289. VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0);
  290. VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols);
  291. RealScalar a = m1.abs().mean();
  292. VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count());
  293. typedef Array<Index, Dynamic, 1> ArrayOfIndices;
  294. // TODO allows colwise/rowwise for array
  295. VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
  296. VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
  297. }
  298. template<typename ArrayType> void array_real(const ArrayType& m)
  299. {
  300. using std::abs;
  301. using std::sqrt;
  302. typedef typename ArrayType::Scalar Scalar;
  303. typedef typename NumTraits<Scalar>::Real RealScalar;
  304. Index rows = m.rows();
  305. Index cols = m.cols();
  306. ArrayType m1 = ArrayType::Random(rows, cols),
  307. m2 = ArrayType::Random(rows, cols),
  308. m3(rows, cols),
  309. m4 = m1;
  310. m4 = (m4.abs()==Scalar(0)).select(Scalar(1),m4);
  311. Scalar s1 = internal::random<Scalar>();
  312. // these tests are mostly to check possible compilation issues with free-functions.
  313. VERIFY_IS_APPROX(m1.sin(), sin(m1));
  314. VERIFY_IS_APPROX(m1.cos(), cos(m1));
  315. VERIFY_IS_APPROX(m1.tan(), tan(m1));
  316. VERIFY_IS_APPROX(m1.asin(), asin(m1));
  317. VERIFY_IS_APPROX(m1.acos(), acos(m1));
  318. VERIFY_IS_APPROX(m1.atan(), atan(m1));
  319. VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
  320. VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
  321. VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
  322. #if EIGEN_HAS_CXX11_MATH
  323. VERIFY_IS_APPROX(m1.tanh().atanh(), atanh(tanh(m1)));
  324. VERIFY_IS_APPROX(m1.sinh().asinh(), asinh(sinh(m1)));
  325. VERIFY_IS_APPROX(m1.cosh().acosh(), acosh(cosh(m1)));
  326. #endif
  327. VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
  328. VERIFY_IS_APPROX(m1.arg(), arg(m1));
  329. VERIFY_IS_APPROX(m1.round(), round(m1));
  330. VERIFY_IS_APPROX(m1.rint(), rint(m1));
  331. VERIFY_IS_APPROX(m1.floor(), floor(m1));
  332. VERIFY_IS_APPROX(m1.ceil(), ceil(m1));
  333. VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
  334. VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
  335. VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
  336. VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
  337. VERIFY_IS_APPROX(m1.abs(), abs(m1));
  338. VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
  339. VERIFY_IS_APPROX(m1.square(), square(m1));
  340. VERIFY_IS_APPROX(m1.cube(), cube(m1));
  341. VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
  342. VERIFY_IS_APPROX(m1.sign(), sign(m1));
  343. VERIFY((m1.sqrt().sign().isNaN() == (Eigen::isnan)(sign(sqrt(m1)))).all());
  344. // avoid inf and NaNs so verification doesn't fail
  345. m3 = m4.abs();
  346. VERIFY_IS_APPROX(m3.sqrt(), sqrt(abs(m3)));
  347. VERIFY_IS_APPROX(m3.rsqrt(), Scalar(1)/sqrt(abs(m3)));
  348. VERIFY_IS_APPROX(rsqrt(m3), Scalar(1)/sqrt(abs(m3)));
  349. VERIFY_IS_APPROX(m3.log(), log(m3));
  350. VERIFY_IS_APPROX(m3.log1p(), log1p(m3));
  351. VERIFY_IS_APPROX(m3.log10(), log10(m3));
  352. VERIFY_IS_APPROX(m3.log2(), log2(m3));
  353. VERIFY((!(m1>m2) == (m1<=m2)).all());
  354. VERIFY_IS_APPROX(sin(m1.asin()), m1);
  355. VERIFY_IS_APPROX(cos(m1.acos()), m1);
  356. VERIFY_IS_APPROX(tan(m1.atan()), m1);
  357. VERIFY_IS_APPROX(sinh(m1), Scalar(0.5)*(exp(m1)-exp(-m1)));
  358. VERIFY_IS_APPROX(cosh(m1), Scalar(0.5)*(exp(m1)+exp(-m1)));
  359. VERIFY_IS_APPROX(tanh(m1), (Scalar(0.5)*(exp(m1)-exp(-m1)))/(Scalar(0.5)*(exp(m1)+exp(-m1))));
  360. VERIFY_IS_APPROX(logistic(m1), (Scalar(1)/(Scalar(1)+exp(-m1))));
  361. VERIFY_IS_APPROX(arg(m1), ((m1<Scalar(0)).template cast<Scalar>())*Scalar(std::acos(Scalar(-1))));
  362. VERIFY((round(m1) <= ceil(m1) && round(m1) >= floor(m1)).all());
  363. VERIFY((rint(m1) <= ceil(m1) && rint(m1) >= floor(m1)).all());
  364. VERIFY(((ceil(m1) - round(m1)) <= Scalar(0.5) || (round(m1) - floor(m1)) <= Scalar(0.5)).all());
  365. VERIFY(((ceil(m1) - round(m1)) <= Scalar(1.0) && (round(m1) - floor(m1)) <= Scalar(1.0)).all());
  366. VERIFY(((ceil(m1) - rint(m1)) <= Scalar(0.5) || (rint(m1) - floor(m1)) <= Scalar(0.5)).all());
  367. VERIFY(((ceil(m1) - rint(m1)) <= Scalar(1.0) && (rint(m1) - floor(m1)) <= Scalar(1.0)).all());
  368. VERIFY((Eigen::isnan)((m1*Scalar(0))/Scalar(0)).all());
  369. VERIFY((Eigen::isinf)(m4/Scalar(0)).all());
  370. VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*Scalar(0)/Scalar(0))) && (!(Eigen::isfinite)(m4/Scalar(0)))).all());
  371. VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
  372. VERIFY((abs(m1) == m1 || abs(m1) == -m1).all());
  373. VERIFY_IS_APPROX(m3, sqrt(abs2(m3)));
  374. VERIFY_IS_APPROX(m1.absolute_difference(m2), (m1 > m2).select(m1 - m2, m2 - m1));
  375. VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
  376. VERIFY_IS_APPROX( m1*m1.sign(),m1.abs());
  377. VERIFY_IS_APPROX(m1.sign() * m1.abs(), m1);
  378. VERIFY_IS_APPROX(numext::abs2(numext::real(m1)) + numext::abs2(numext::imag(m1)), numext::abs2(m1));
  379. VERIFY_IS_APPROX(numext::abs2(Eigen::real(m1)) + numext::abs2(Eigen::imag(m1)), numext::abs2(m1));
  380. if(!NumTraits<Scalar>::IsComplex)
  381. VERIFY_IS_APPROX(numext::real(m1), m1);
  382. // shift argument of logarithm so that it is not zero
  383. Scalar smallNumber = NumTraits<Scalar>::dummy_precision();
  384. VERIFY_IS_APPROX((m3 + smallNumber).log() , log(abs(m3) + smallNumber));
  385. VERIFY_IS_APPROX((m3 + smallNumber + Scalar(1)).log() , log1p(abs(m3) + smallNumber));
  386. VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
  387. VERIFY_IS_APPROX(m1.exp(), exp(m1));
  388. VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
  389. VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
  390. VERIFY_IS_APPROX((m3 + smallNumber).exp() - Scalar(1), expm1(abs(m3) + smallNumber));
  391. VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
  392. VERIFY_IS_APPROX(pow(m3,RealScalar(0.5)), m3.sqrt());
  393. VERIFY_IS_APPROX(m3.pow(RealScalar(-0.5)), m3.rsqrt());
  394. VERIFY_IS_APPROX(pow(m3,RealScalar(-0.5)), m3.rsqrt());
  395. // Avoid inf and NaN.
  396. m3 = (m1.square()<NumTraits<Scalar>::epsilon()).select(Scalar(1),m3);
  397. VERIFY_IS_APPROX(m3.pow(RealScalar(-2)), m3.square().inverse());
  398. pow_test<Scalar>();
  399. VERIFY_IS_APPROX(log10(m3), log(m3)/numext::log(Scalar(10)));
  400. VERIFY_IS_APPROX(log2(m3), log(m3)/numext::log(Scalar(2)));
  401. // scalar by array division
  402. const RealScalar tiny = sqrt(std::numeric_limits<RealScalar>::epsilon());
  403. s1 += Scalar(tiny);
  404. m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
  405. VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
  406. // check inplace transpose
  407. m3 = m1;
  408. m3.transposeInPlace();
  409. VERIFY_IS_APPROX(m3, m1.transpose());
  410. m3.transposeInPlace();
  411. VERIFY_IS_APPROX(m3, m1);
  412. }
  413. template<typename ArrayType> void array_complex(const ArrayType& m)
  414. {
  415. typedef typename ArrayType::Scalar Scalar;
  416. typedef typename NumTraits<Scalar>::Real RealScalar;
  417. Index rows = m.rows();
  418. Index cols = m.cols();
  419. ArrayType m1 = ArrayType::Random(rows, cols),
  420. m2(rows, cols),
  421. m4 = m1;
  422. m4.real() = (m4.real().abs()==RealScalar(0)).select(RealScalar(1),m4.real());
  423. m4.imag() = (m4.imag().abs()==RealScalar(0)).select(RealScalar(1),m4.imag());
  424. Array<RealScalar, -1, -1> m3(rows, cols);
  425. for (Index i = 0; i < m.rows(); ++i)
  426. for (Index j = 0; j < m.cols(); ++j)
  427. m2(i,j) = sqrt(m1(i,j));
  428. // these tests are mostly to check possible compilation issues with free-functions.
  429. VERIFY_IS_APPROX(m1.sin(), sin(m1));
  430. VERIFY_IS_APPROX(m1.cos(), cos(m1));
  431. VERIFY_IS_APPROX(m1.tan(), tan(m1));
  432. VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
  433. VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
  434. VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
  435. VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
  436. VERIFY_IS_APPROX(m1.arg(), arg(m1));
  437. VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
  438. VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
  439. VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
  440. VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
  441. VERIFY_IS_APPROX(m1.log(), log(m1));
  442. VERIFY_IS_APPROX(m1.log10(), log10(m1));
  443. VERIFY_IS_APPROX(m1.log2(), log2(m1));
  444. VERIFY_IS_APPROX(m1.abs(), abs(m1));
  445. VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
  446. VERIFY_IS_APPROX(m1.sqrt(), sqrt(m1));
  447. VERIFY_IS_APPROX(m1.square(), square(m1));
  448. VERIFY_IS_APPROX(m1.cube(), cube(m1));
  449. VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
  450. VERIFY_IS_APPROX(m1.sign(), sign(m1));
  451. VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
  452. VERIFY_IS_APPROX(m1.exp(), exp(m1));
  453. VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
  454. VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
  455. VERIFY_IS_APPROX(expm1(m1), exp(m1) - 1.);
  456. // Check for larger magnitude complex numbers that expm1 matches exp - 1.
  457. VERIFY_IS_APPROX(expm1(10. * m1), exp(10. * m1) - 1.);
  458. VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
  459. VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
  460. VERIFY_IS_APPROX(tanh(m1), (0.5*(exp(m1)-exp(-m1)))/(0.5*(exp(m1)+exp(-m1))));
  461. VERIFY_IS_APPROX(logistic(m1), (1.0/(1.0 + exp(-m1))));
  462. for (Index i = 0; i < m.rows(); ++i)
  463. for (Index j = 0; j < m.cols(); ++j)
  464. m3(i,j) = std::atan2(m1(i,j).imag(), m1(i,j).real());
  465. VERIFY_IS_APPROX(arg(m1), m3);
  466. std::complex<RealScalar> zero(0.0,0.0);
  467. VERIFY((Eigen::isnan)(m1*zero/zero).all());
  468. #if EIGEN_COMP_MSVC
  469. // msvc complex division is not robust
  470. VERIFY((Eigen::isinf)(m4/RealScalar(0)).all());
  471. #else
  472. #if EIGEN_COMP_CLANG
  473. // clang's complex division is notoriously broken too
  474. if((numext::isinf)(m4(0,0)/RealScalar(0))) {
  475. #endif
  476. VERIFY((Eigen::isinf)(m4/zero).all());
  477. #if EIGEN_COMP_CLANG
  478. }
  479. else
  480. {
  481. VERIFY((Eigen::isinf)(m4.real()/zero.real()).all());
  482. }
  483. #endif
  484. #endif // MSVC
  485. VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*zero/zero)) && (!(Eigen::isfinite)(m1/zero))).all());
  486. VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
  487. VERIFY_IS_APPROX(conj(m1.conjugate()), m1);
  488. VERIFY_IS_APPROX(abs(m1), sqrt(square(m1.real())+square(m1.imag())));
  489. VERIFY_IS_APPROX(abs(m1), sqrt(abs2(m1)));
  490. VERIFY_IS_APPROX(log10(m1), log(m1)/log(10));
  491. VERIFY_IS_APPROX(log2(m1), log(m1)/log(2));
  492. VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
  493. VERIFY_IS_APPROX( m1.sign() * m1.abs(), m1);
  494. // scalar by array division
  495. Scalar s1 = internal::random<Scalar>();
  496. const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon());
  497. s1 += Scalar(tiny);
  498. m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
  499. VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
  500. // check inplace transpose
  501. m2 = m1;
  502. m2.transposeInPlace();
  503. VERIFY_IS_APPROX(m2, m1.transpose());
  504. m2.transposeInPlace();
  505. VERIFY_IS_APPROX(m2, m1);
  506. // Check vectorized inplace transpose.
  507. ArrayType m5 = ArrayType::Random(131, 131);
  508. ArrayType m6 = m5;
  509. m6.transposeInPlace();
  510. VERIFY_IS_APPROX(m6, m5.transpose());
  511. }
  512. template<typename ArrayType> void min_max(const ArrayType& m)
  513. {
  514. typedef typename ArrayType::Scalar Scalar;
  515. Index rows = m.rows();
  516. Index cols = m.cols();
  517. ArrayType m1 = ArrayType::Random(rows, cols);
  518. // min/max with array
  519. Scalar maxM1 = m1.maxCoeff();
  520. Scalar minM1 = m1.minCoeff();
  521. VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1)));
  522. VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1)));
  523. VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1)));
  524. VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1)));
  525. // min/max with scalar input
  526. VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1));
  527. VERIFY_IS_APPROX(m1, (m1.min)( maxM1));
  528. VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1));
  529. VERIFY_IS_APPROX(m1, (m1.max)( minM1));
  530. // min/max with various NaN propagation options.
  531. if (m1.size() > 1 && !NumTraits<Scalar>::IsInteger) {
  532. m1(0,0) = NumTraits<Scalar>::quiet_NaN();
  533. maxM1 = m1.template maxCoeff<PropagateNaN>();
  534. minM1 = m1.template minCoeff<PropagateNaN>();
  535. VERIFY((numext::isnan)(maxM1));
  536. VERIFY((numext::isnan)(minM1));
  537. maxM1 = m1.template maxCoeff<PropagateNumbers>();
  538. minM1 = m1.template minCoeff<PropagateNumbers>();
  539. VERIFY(!(numext::isnan)(maxM1));
  540. VERIFY(!(numext::isnan)(minM1));
  541. }
  542. }
  543. template<int N>
  544. struct shift_left {
  545. template<typename Scalar>
  546. Scalar operator()(const Scalar& v) const {
  547. return v << N;
  548. }
  549. };
  550. template<int N>
  551. struct arithmetic_shift_right {
  552. template<typename Scalar>
  553. Scalar operator()(const Scalar& v) const {
  554. return v >> N;
  555. }
  556. };
  557. template<typename ArrayType> void array_integer(const ArrayType& m)
  558. {
  559. Index rows = m.rows();
  560. Index cols = m.cols();
  561. ArrayType m1 = ArrayType::Random(rows, cols),
  562. m2(rows, cols);
  563. m2 = m1.template shiftLeft<2>();
  564. VERIFY( (m2 == m1.unaryExpr(shift_left<2>())).all() );
  565. m2 = m1.template shiftLeft<9>();
  566. VERIFY( (m2 == m1.unaryExpr(shift_left<9>())).all() );
  567. m2 = m1.template shiftRight<2>();
  568. VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<2>())).all() );
  569. m2 = m1.template shiftRight<9>();
  570. VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<9>())).all() );
  571. }
  572. EIGEN_DECLARE_TEST(array_cwise)
  573. {
  574. for(int i = 0; i < g_repeat; i++) {
  575. CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
  576. CALL_SUBTEST_2( array(Array22f()) );
  577. CALL_SUBTEST_3( array(Array44d()) );
  578. CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  579. CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  580. CALL_SUBTEST_6( array(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  581. CALL_SUBTEST_6( array(Array<Index,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  582. CALL_SUBTEST_6( array_integer(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  583. CALL_SUBTEST_6( array_integer(Array<Index,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  584. }
  585. for(int i = 0; i < g_repeat; i++) {
  586. CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
  587. CALL_SUBTEST_2( comparisons(Array22f()) );
  588. CALL_SUBTEST_3( comparisons(Array44d()) );
  589. CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  590. CALL_SUBTEST_6( comparisons(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  591. }
  592. for(int i = 0; i < g_repeat; i++) {
  593. CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) );
  594. CALL_SUBTEST_2( min_max(Array22f()) );
  595. CALL_SUBTEST_3( min_max(Array44d()) );
  596. CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  597. CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  598. }
  599. for(int i = 0; i < g_repeat; i++) {
  600. CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
  601. CALL_SUBTEST_2( array_real(Array22f()) );
  602. CALL_SUBTEST_3( array_real(Array44d()) );
  603. CALL_SUBTEST_5( array_real(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  604. CALL_SUBTEST_7( array_real(Array<Eigen::half, 32, 32>()) );
  605. CALL_SUBTEST_8( array_real(Array<Eigen::bfloat16, 32, 32>()) );
  606. }
  607. for(int i = 0; i < g_repeat; i++) {
  608. CALL_SUBTEST_4( array_complex(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
  609. }
  610. VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
  611. VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
  612. VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
  613. typedef CwiseUnaryOp<internal::scalar_abs_op<double>, ArrayXd > Xpr;
  614. VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
  615. ArrayBase<Xpr>
  616. >::value));
  617. }