blob: a7e0ff48cae387a2dad3fc65e633ecb1e5e2320e [file] [log] [blame]
Gael Guennebaud8cef5412008-06-21 17:28:07 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaud8cef5412008-06-21 17:28:07 +00003//
Gael Guennebaud28e64b02010-06-24 23:21:58 +02004// Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr>
Gael Guennebaud8cef5412008-06-21 17:28:07 +00005//
Benoit Jacob69124cf2012-07-13 14:42:47 -04006// 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/.
Gael Guennebaud19bc4182013-02-25 01:17:44 +01009
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000010#include <vector>
Gael Guennebaud8cef5412008-06-21 17:28:07 +000011#include "main.h"
Gael Guennebaud8cef5412008-06-21 17:28:07 +000012
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000013template <typename Scalar>
14std::vector<Scalar> special_values() {
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000015 const Scalar zero = Scalar(0);
Antonio Sanchez8dfe1022021-03-16 20:12:46 -070016 const Scalar eps = Eigen::NumTraits<Scalar>::epsilon();
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000017 const Scalar one = Scalar(1);
Rasmus Munk Larsen6e3b7952021-02-05 16:58:49 -080018 const Scalar two = Scalar(2);
19 const Scalar three = Scalar(3);
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000020 const Scalar sqrt_half = Scalar(std::sqrt(0.5));
21 const Scalar sqrt2 = Scalar(std::sqrt(2));
Antonio Sanchez8dfe1022021-03-16 20:12:46 -070022 const Scalar inf = Eigen::NumTraits<Scalar>::infinity();
23 const Scalar nan = Eigen::NumTraits<Scalar>::quiet_NaN();
Rasmus Munk Larsenbe0574e2021-02-17 02:50:32 +000024 const Scalar denorm_min = std::numeric_limits<Scalar>::denorm_min();
Rasmus Munk Larsen6e3b7952021-02-05 16:58:49 -080025 const Scalar min = (std::numeric_limits<Scalar>::min)();
26 const Scalar max = (std::numeric_limits<Scalar>::max)();
Antonio Sanchez8dfe1022021-03-16 20:12:46 -070027 const Scalar max_exp = (static_cast<Scalar>(int(Eigen::NumTraits<Scalar>::max_exponent())) * Scalar(EIGEN_LN2)) / eps;
Rasmus Munk Larsenbe0574e2021-02-17 02:50:32 +000028
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000029 return {zero, denorm_min, min, eps, sqrt_half, one, sqrt2, two, three, max_exp, max, inf, nan};
30}
31
32template<typename Scalar>
33void special_value_pairs(Array<Scalar, Dynamic, Dynamic>& x,
34 Array<Scalar, Dynamic, Dynamic>& y) {
35 std::vector<Scalar> abs_vals = special_values<Scalar>();
36 const int abs_cases = abs_vals.size();
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000037 const int num_cases = 2*abs_cases * 2*abs_cases;
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000038 // ensure both vectorized and non-vectorized paths taken
39 const int num_repeats = 2 * internal::packet_traits<Scalar>::size + 1;
40 x.resize(num_repeats, num_cases);
41 y.resize(num_repeats, num_cases);
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000042 int count = 0;
43 for (int i = 0; i < abs_cases; ++i) {
44 const Scalar abs_x = abs_vals[i];
45 for (int sign_x = 0; sign_x < 2; ++sign_x) {
46 Scalar x_case = sign_x == 0 ? -abs_x : abs_x;
47 for (int j = 0; j < abs_cases; ++j) {
48 const Scalar abs_y = abs_vals[j];
49 for (int sign_y = 0; sign_y < 2; ++sign_y) {
50 Scalar y_case = sign_y == 0 ? -abs_y : abs_y;
51 for (int repeat = 0; repeat < num_repeats; ++repeat) {
52 x(repeat, count) = x_case;
53 y(repeat, count) = y_case;
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000054 }
55 ++count;
56 }
57 }
58 }
59 }
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000060}
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000061
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000062template <typename Scalar, typename Fn, typename RefFn>
63void binary_op_test(std::string name, Fn fun, RefFn ref) {
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000064 const Scalar tol = test_precision<Scalar>();
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000065 Array<Scalar, Dynamic, Dynamic> x;
66 Array<Scalar, Dynamic, Dynamic> y;
67 special_value_pairs(x, y);
68
69 Array<Scalar, Dynamic, Dynamic> actual = fun(x, y);
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000070 bool all_pass = true;
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000071 for (int i = 0; i < x.rows(); ++i) {
72 for (int j = 0; j < x.cols(); ++j) {
73 Scalar e = static_cast<Scalar>(ref(x(i,j), y(i,j)));
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000074 Scalar a = actual(i, j);
Rasmus Munk Larsen7a87ed12022-08-08 18:48:36 +000075 bool success = (a==e) || ((numext::isfinite)(e) && internal::isApprox(a, e, tol)) || ((numext::isnan)(a) && (numext::isnan)(e));
76 all_pass &= success;
77 if (!success) {
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000078 std::cout << name << "(" << x(i,j) << "," << y(i,j) << ") = " << a << " != " << e << std::endl;
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000079 }
80 }
81 }
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000082 VERIFY(all_pass);
83}
Charles Schlosser76a669f2022-08-16 21:32:36 +000084
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000085template <typename Scalar>
86void binary_ops_test() {
87 binary_op_test<Scalar>("pow",
88 [](auto x, auto y) { return Eigen::pow(x, y); },
89 [](auto x, auto y) { return std::pow(x, y); });
90 binary_op_test<Scalar>("atan2",
91 [](auto x, auto y) { return Eigen::atan2(x, y); },
92 [](auto x, auto y) { return std::atan2(x, y); });
93}
Charles Schlosser76a669f2022-08-16 21:32:36 +000094
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +000095template <typename Scalar>
96void pow_scalar_exponent_test() {
97 using Int_t = typename internal::make_integer<Scalar>::type;
98 const Scalar tol = test_precision<Scalar>();
99
100 std::vector<Scalar> abs_vals = special_values<Scalar>();
101 const int num_vals = abs_vals.size();
102 Map<Array<Scalar, Dynamic, 1>> bases(abs_vals.data(), num_vals);
103
104 bool all_pass = true;
105 for (Scalar abs_exponent : abs_vals) {
106 for (Scalar exponent : {-abs_exponent, abs_exponent}) {
107 // test integer exponent code path
108 bool exponent_is_integer = (numext::isfinite)(exponent) && (numext::round(exponent) == exponent) &&
109 (numext::abs(exponent) < static_cast<Scalar>(NumTraits<Int_t>::highest()));
110 if (exponent_is_integer) {
111 Int_t exponent_as_int = static_cast<Int_t>(exponent);
112 Array<Scalar, Dynamic, 1> eigenPow = bases.pow(exponent_as_int);
113 for (int j = 0; j < num_vals; j++) {
Charles Schlosser76a669f2022-08-16 21:32:36 +0000114 Scalar e = static_cast<Scalar>(std::pow(bases(j), exponent));
115 Scalar a = eigenPow(j);
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +0000116 bool success = (a == e) || ((numext::isfinite)(e) && internal::isApprox(a, e, tol)) ||
117 ((numext::isnan)(a) && (numext::isnan)(e));
Charles Schlosser76a669f2022-08-16 21:32:36 +0000118 all_pass &= success;
119 if (!success) {
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +0000120 std::cout << "pow(" << bases(j) << "," << exponent << ") = " << a << " != " << e << std::endl;
Charles Schlosser76a669f2022-08-16 21:32:36 +0000121 }
122 }
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +0000123 } else {
124 // test floating point exponent code path
125 Array<Scalar, Dynamic, 1> eigenPow = bases.pow(exponent);
126 for (int j = 0; j < num_vals; j++) {
127 Scalar e = static_cast<Scalar>(std::pow(bases(j), exponent));
128 Scalar a = eigenPow(j);
129 bool success = (a == e) || ((numext::isfinite)(e) && internal::isApprox(a, e, tol)) ||
130 ((numext::isnan)(a) && (numext::isnan)(e));
131 all_pass &= success;
132 if (!success) {
133 std::cout << "pow(" << bases(j) << "," << exponent << ") = " << a << " != " << e << std::endl;
Charles Schlosser76a669f2022-08-16 21:32:36 +0000134 }
135 }
136 }
137 }
138 }
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +0000139 VERIFY(all_pass);
140}
141
Charles Schlossere5af9f82022-08-29 19:23:54 +0000142template <typename Scalar, typename ScalarExponent>
143Scalar calc_overflow_threshold(const ScalarExponent exponent) {
144 EIGEN_USING_STD(exp2);
Antonio Sánchez69f50e32022-09-06 19:53:29 +0000145 EIGEN_USING_STD(log2);
Charles Schlossere5af9f82022-08-29 19:23:54 +0000146 EIGEN_STATIC_ASSERT((NumTraits<Scalar>::digits() < 2 * NumTraits<double>::digits()), BASE_TYPE_IS_TOO_BIG);
147
148 if (exponent < 2)
149 return NumTraits<Scalar>::highest();
150 else {
Antonio Sánchez69f50e32022-09-06 19:53:29 +0000151 // base^e <= highest ==> base <= 2^(log2(highest)/e)
Antonio Sánchez3c37dd22022-09-08 17:40:45 +0000152 // For floating-point types, consider the bound for integer values that can be reproduced exactly = 2 ^ digits
153 double highest_bits = numext::mini(static_cast<double>(NumTraits<Scalar>::digits()),
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000154 static_cast<double>(log2(NumTraits<Scalar>::highest())));
Antonio Sánchez69f50e32022-09-06 19:53:29 +0000155 return static_cast<Scalar>(
Antonio Sánchez3c37dd22022-09-08 17:40:45 +0000156 numext::floor(exp2(highest_bits / static_cast<double>(exponent))));
Charles Schlossere5af9f82022-08-29 19:23:54 +0000157 }
158}
159
Rasmus Munk Larsendceb7792022-09-12 15:51:27 -0700160template <typename Base, typename Exponent, bool ExpIsInteger = NumTraits<Exponent>::IsInteger>
161struct ref_pow {
162 static Base run(Base base, Exponent exponent) {
163 EIGEN_USING_STD(pow);
164 return pow(base, static_cast<Base>(exponent));
165 }
166};
167
168template <typename Base, typename Exponent>
169struct ref_pow<Base, Exponent, true> {
170 static Base run(Base base, Exponent exponent) {
171 EIGEN_USING_STD(pow);
172 return pow(base, exponent);
173 }
174};
175
Charles Schlossere5af9f82022-08-29 19:23:54 +0000176template <typename Base, typename Exponent>
177void test_exponent(Exponent exponent) {
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000178 const Base max_abs_bases = static_cast<Base>(10000);
179 // avoid integer overflow in Base type
180 Base threshold = calc_overflow_threshold<Base, Exponent>(numext::abs(exponent));
181 // avoid numbers that can't be verified with std::pow
182 double double_threshold = calc_overflow_threshold<double, Exponent>(numext::abs(exponent));
183 // use the lesser of these two thresholds
184 Base testing_threshold =
185 static_cast<double>(threshold) < double_threshold ? threshold : static_cast<Base>(double_threshold);
186 // test both vectorized and non-vectorized code paths
187 const Index array_size = 2 * internal::packet_traits<Base>::size + 1;
188
189 Base max_base = numext::mini(testing_threshold, max_abs_bases);
190 Base min_base = NumTraits<Base>::IsSigned ? -max_base : Base(0);
191
192 ArrayX<Base> x(array_size), y(array_size);
193 bool all_pass = true;
194 for (Base base = min_base; base <= max_base; base++) {
195 if (exponent < 0 && base == 0) continue;
196 x.setConstant(base);
197 y = x.pow(exponent);
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000198 for (Base a : y) {
Rasmus Munk Larsendceb7792022-09-12 15:51:27 -0700199 Base e = ref_pow<Base, Exponent>::run(base, exponent);
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000200 bool pass = (a == e);
201 if (!NumTraits<Base>::IsInteger) {
202 pass = pass || (((numext::isfinite)(e) && internal::isApprox(a, e)) ||
203 ((numext::isnan)(a) && (numext::isnan)(e)));
204 }
205 all_pass &= pass;
206 if (!pass) {
207 std::cout << "pow(" << base << "," << exponent << ") = " << a << " != " << e << std::endl;
208 }
Charles Schlossere5af9f82022-08-29 19:23:54 +0000209 }
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000210 }
211 VERIFY(all_pass);
Charles Schlossere5af9f82022-08-29 19:23:54 +0000212}
Charles Schlossere5af9f82022-08-29 19:23:54 +0000213
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000214template <typename Base, typename Exponent>
215void unary_pow_test() {
216 Exponent max_exponent = static_cast<Exponent>(NumTraits<Base>::digits());
217 Exponent min_exponent = static_cast<Exponent>(NumTraits<Exponent>::IsSigned ? -max_exponent : 0);
218
219 for (Exponent exponent = min_exponent; exponent < max_exponent; ++exponent) {
220 test_exponent<Base, Exponent>(exponent);
221 }
222};
223
224void mixed_pow_test() {
225 // The following cases will test promoting a smaller exponent type
226 // to a wider base type.
227 unary_pow_test<double, int>();
228 unary_pow_test<double, float>();
229 unary_pow_test<float, half>();
230 unary_pow_test<double, half>();
231 unary_pow_test<float, bfloat16>();
232 unary_pow_test<double, bfloat16>();
233
234 // Although in the following cases the exponent cannot be represented exactly
235 // in the base type, we do not perform a conversion, but implement
236 // the operation using repeated squaring.
237 unary_pow_test<float, int>();
238 unary_pow_test<double, long long>();
239
240 // The following cases will test promoting a wider exponent type
241 // to a narrower base type. This should compile but generate a
242 // deprecation warning:
243 unary_pow_test<float, double>();
244}
245
246void int_pow_test() {
247 unary_pow_test<int, int>();
248 unary_pow_test<unsigned int, unsigned int>();
249 unary_pow_test<long long, long long>();
250 unary_pow_test<unsigned long long, unsigned long long>();
251
252 // Although in the following cases the exponent cannot be represented exactly
253 // in the base type, we do not perform a conversion, but implement the
254 // operation using repeated squaring.
255 unary_pow_test<long long, int>();
256 unary_pow_test<int, unsigned int>();
257 unary_pow_test<unsigned int, int>();
258 unary_pow_test<long long, unsigned long long>();
259 unary_pow_test<unsigned long long, long long>();
260 unary_pow_test<long long, int>();
Charles Schlossere5af9f82022-08-29 19:23:54 +0000261}
262
Benoit Jacobe2775862010-04-28 18:51:38 -0400263template<typename ArrayType> void array(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000264{
Benoit Jacobe2775862010-04-28 18:51:38 -0400265 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200266 typedef typename ArrayType::RealScalar RealScalar;
Benoit Jacobe2775862010-04-28 18:51:38 -0400267 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
268 typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000269
Hauke Heibelf1679c72010-06-20 17:37:56 +0200270 Index rows = m.rows();
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800271 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000272
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000273 ArrayType m1 = ArrayType::Random(rows, cols);
274 if (NumTraits<RealScalar>::IsInteger && NumTraits<RealScalar>::IsSigned
275 && !NumTraits<Scalar>::IsComplex) {
276 // Here we cap the size of the values in m1 such that pow(3)/cube()
277 // doesn't overflow and result in undefined behavior. Notice that because
278 // pow(int, int) promotes its inputs and output to double (according to
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000279 // the C++ standard), we have to make sure that the result fits in 53 bits
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000280 // for int64,
281 RealScalar max_val =
282 numext::mini(RealScalar(std::cbrt(NumTraits<RealScalar>::highest())),
283 RealScalar(std::cbrt(1LL << 53)))/2;
284 m1.array() = (m1.abs().array() <= max_val).select(m1, Scalar(max_val));
285 }
286 ArrayType m2 = ArrayType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000287 m3(rows, cols);
Christoph Hertzberg3be9f5c2015-04-16 13:25:20 +0200288 ArrayType m4 = m1; // copy constructor
289 VERIFY_IS_APPROX(m1, m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000290
Gael Guennebaud627595a2009-06-10 11:20:30 +0200291 ColVectorType cv1 = ColVectorType::Random(rows);
292 RowVectorType rv1 = RowVectorType::Random(cols);
293
Benoit Jacob47160402010-10-25 10:15:22 -0400294 Scalar s1 = internal::random<Scalar>(),
Hauke Heibel8cb3e362012-03-02 16:27:27 +0100295 s2 = internal::random<Scalar>();
296
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000297 // scalar addition
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100298 VERIFY_IS_APPROX(m1 + s1, s1 + m1);
Benoit Jacobe2775862010-04-28 18:51:38 -0400299 VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100300 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
Benoit Jacobe2775862010-04-28 18:51:38 -0400301 VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
302 VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
303 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000304 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100305 m3 += s2;
306 VERIFY_IS_APPROX(m3, m1 + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000307 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100308 m3 -= s1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800309 VERIFY_IS_APPROX(m3, m1 - s1);
310
Hauke Heibelf578dc72010-12-16 17:34:13 +0100311 // scalar operators via Maps
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000312 m3 = m1; m4 = m1;
313 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
314 VERIFY_IS_APPROX(m4, m3 - m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800315
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000316 m3 = m1; m4 = m1;
317 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
318 VERIFY_IS_APPROX(m4, m3 + m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800319
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000320 m3 = m1; m4 = m1;
321 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
322 VERIFY_IS_APPROX(m4, m3 * m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800323
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000324 m3 = m1; m4 = m1;
Hauke Heibelf578dc72010-12-16 17:34:13 +0100325 m2 = ArrayType::Random(rows,cols);
326 m2 = (m2==0).select(1,m2);
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000327 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
328 VERIFY_IS_APPROX(m4, m3 / m2);
Gael Guennebaud05ad0832008-07-19 13:03:23 +0000329
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000330 // reductions
Gael Guennebaud42af5872013-02-23 22:58:14 +0100331 VERIFY_IS_APPROX(m1.abs().colwise().sum().sum(), m1.abs().sum());
332 VERIFY_IS_APPROX(m1.abs().rowwise().sum().sum(), m1.abs().sum());
333 using std::abs;
334 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.colwise().sum().sum() - m1.sum()), m1.abs().sum());
335 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.rowwise().sum().sum() - m1.sum()), m1.abs().sum());
Gael Guennebaud28e139a2013-02-23 23:06:45 +0100336 if (!internal::isMuchSmallerThan(abs(m1.sum() - (m1+m2).sum()), m1.abs().sum(), test_precision<Scalar>()))
Hauke Heibel83e3c452010-12-16 18:53:02 +0100337 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud66e99ab2016-06-06 15:11:41 +0200338 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +0200339
340 // vector-wise ops
341 m3 = m1;
342 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
343 m3 = m1;
344 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
345 m3 = m1;
346 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
347 m3 = m1;
348 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800349
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200350 // Conversion from scalar
351 VERIFY_IS_APPROX((m3 = s1), ArrayType::Constant(rows,cols,s1));
352 VERIFY_IS_APPROX((m3 = 1), ArrayType::Constant(rows,cols,1));
353 VERIFY_IS_APPROX((m3.topLeftCorner(rows,cols) = 1), ArrayType::Constant(rows,cols,1));
354 typedef Array<Scalar,
355 ArrayType::RowsAtCompileTime==Dynamic?2:ArrayType::RowsAtCompileTime,
356 ArrayType::ColsAtCompileTime==Dynamic?2:ArrayType::ColsAtCompileTime,
357 ArrayType::Options> FixedArrayType;
Gael Guennebaud543529d2019-01-22 15:30:50 +0100358 {
359 FixedArrayType f1(s1);
360 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
361 FixedArrayType f2(numext::real(s1));
362 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
363 FixedArrayType f3((int)100*numext::real(s1));
364 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
365 f1.setRandom();
366 FixedArrayType f4(f1.data());
367 VERIFY_IS_APPROX(f4, f1);
368 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100369 {
370 FixedArrayType f1{s1};
371 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
372 FixedArrayType f2{numext::real(s1)};
373 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
374 FixedArrayType f3{(int)100*numext::real(s1)};
375 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
376 f1.setRandom();
377 FixedArrayType f4{f1.data()};
378 VERIFY_IS_APPROX(f4, f1);
379 }
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800380
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200381 // pow
382 VERIFY_IS_APPROX(m1.pow(2), m1.square());
383 VERIFY_IS_APPROX(pow(m1,2), m1.square());
384 VERIFY_IS_APPROX(m1.pow(3), m1.cube());
385 VERIFY_IS_APPROX(pow(m1,3), m1.cube());
386 VERIFY_IS_APPROX((-m1).pow(3), -m1.cube());
387 VERIFY_IS_APPROX(pow(2*m1,3), 8*m1.cube());
388 ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2));
389 VERIFY_IS_APPROX(Eigen::pow(m1,exponents), m1.square());
390 VERIFY_IS_APPROX(m1.pow(exponents), m1.square());
391 VERIFY_IS_APPROX(Eigen::pow(2*m1,exponents), 4*m1.square());
392 VERIFY_IS_APPROX((2*m1).pow(exponents), 4*m1.square());
393 VERIFY_IS_APPROX(Eigen::pow(m1,2*exponents), m1.square().square());
394 VERIFY_IS_APPROX(m1.pow(2*exponents), m1.square().square());
Gael Guennebaude61cee72016-07-04 11:49:03 +0200395 VERIFY_IS_APPROX(Eigen::pow(m1(0,0), exponents), ArrayType::Constant(rows,cols,m1(0,0)*m1(0,0)));
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200396
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200397 // Check possible conflicts with 1D ctor
398 typedef Array<Scalar, Dynamic, 1> OneDArrayType;
Gael Guennebaud543529d2019-01-22 15:30:50 +0100399 {
400 OneDArrayType o1(rows);
401 VERIFY(o1.size()==rows);
402 OneDArrayType o2(static_cast<int>(rows));
403 VERIFY(o2.size()==rows);
404 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100405 {
406 OneDArrayType o1{rows};
407 VERIFY(o1.size()==rows);
408 OneDArrayType o4{int(rows)};
409 VERIFY(o4.size()==rows);
410 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100411 // Check possible conflicts with 2D ctor
412 typedef Array<Scalar, Dynamic, Dynamic> TwoDArrayType;
413 typedef Array<Scalar, 2, 1> ArrayType2;
414 {
415 TwoDArrayType o1(rows,cols);
416 VERIFY(o1.rows()==rows);
417 VERIFY(o1.cols()==cols);
418 TwoDArrayType o2(static_cast<int>(rows),static_cast<int>(cols));
419 VERIFY(o2.rows()==rows);
420 VERIFY(o2.cols()==cols);
421
422 ArrayType2 o3(rows,cols);
423 VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
424 ArrayType2 o4(static_cast<int>(rows),static_cast<int>(cols));
425 VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
426 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100427 {
428 TwoDArrayType o1{rows,cols};
429 VERIFY(o1.rows()==rows);
430 VERIFY(o1.cols()==cols);
431 TwoDArrayType o2{int(rows),int(cols)};
432 VERIFY(o2.rows()==rows);
433 VERIFY(o2.cols()==cols);
434
435 ArrayType2 o3{rows,cols};
436 VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
437 ArrayType2 o4{int(rows),int(cols)};
438 VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
439 }
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000440}
441
Benoit Jacobe2775862010-04-28 18:51:38 -0400442template<typename ArrayType> void comparisons(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000443{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100444 using std::abs;
Benoit Jacobe2775862010-04-28 18:51:38 -0400445 typedef typename ArrayType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +0000446 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000447
Hauke Heibelf1679c72010-06-20 17:37:56 +0200448 Index rows = m.rows();
449 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000450
Benoit Jacob47160402010-10-25 10:15:22 -0400451 Index r = internal::random<Index>(0, rows-1),
452 c = internal::random<Index>(0, cols-1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000453
Benoit Jacobe2775862010-04-28 18:51:38 -0400454 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200455 m2 = ArrayType::Random(rows, cols),
456 m3(rows, cols),
457 m4 = m1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800458
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200459 m4 = (m4.abs()==Scalar(0)).select(1,m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000460
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100461 VERIFY(((m1 + Scalar(1)) > m1).all());
462 VERIFY(((m1 - Scalar(1)) < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000463 if (rows*cols>1)
464 {
465 m3 = m1;
466 m3(r,c) += 1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100467 VERIFY(! (m1 < m3).all() );
468 VERIFY(! (m1 > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000469 }
Christoph Hertzberg494fa992015-05-07 17:28:40 +0200470 VERIFY(!(m1 > m2 && m1 < m2).any());
471 VERIFY((m1 <= m2 || m1 >= m2).all());
Gael Guennebaud627595a2009-06-10 11:20:30 +0200472
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200473 // comparisons array to scalar
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100474 VERIFY( (m1 != (m1(r,c)+1) ).any() );
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200475 VERIFY( (m1 > (m1(r,c)-1) ).any() );
476 VERIFY( (m1 < (m1(r,c)+1) ).any() );
477 VERIFY( (m1 == m1(r,c) ).any() );
478
479 // comparisons scalar to array
480 VERIFY( ( (m1(r,c)+1) != m1).any() );
481 VERIFY( ( (m1(r,c)-1) < m1).any() );
482 VERIFY( ( (m1(r,c)+1) > m1).any() );
483 VERIFY( ( m1(r,c) == m1).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200484
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000485 // test Select
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100486 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
487 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
Gael Guennebaud9f795582009-12-16 19:18:40 +0100488 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000489 for (int j=0; j<cols; ++j)
490 for (int i=0; i<rows; ++i)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100491 m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
Benoit Jacobe2775862010-04-28 18:51:38 -0400492 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
493 .select(ArrayType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000494 // shorter versions:
Benoit Jacobe2775862010-04-28 18:51:38 -0400495 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000496 .select(0,m1), m3);
Benoit Jacobe2775862010-04-28 18:51:38 -0400497 VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000498 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000499 // even shorter version:
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100500 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200501
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000502 // count
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100503 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
Benoit Jacobaaaade42010-05-30 16:00:58 -0400504
Gael Guennebaud35c11582011-05-31 22:17:34 +0200505 // and/or
506 VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0);
507 VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols);
508 RealScalar a = m1.abs().mean();
509 VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count());
510
Gael Guennebaud12e1ebb2018-07-12 17:16:40 +0200511 typedef Array<Index, Dynamic, 1> ArrayOfIndices;
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200512
Gael Guennebaud9f795582009-12-16 19:18:40 +0100513 // TODO allows colwise/rowwise for array
Benoit Jacobaaaade42010-05-30 16:00:58 -0400514 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
515 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
Benoit Jacobe8009992008-11-03 22:47:00 +0000516}
517
Benoit Jacobe2775862010-04-28 18:51:38 -0400518template<typename ArrayType> void array_real(const ArrayType& m)
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100519{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100520 using std::abs;
Gael Guennebaud9b9177f2013-07-05 13:35:34 +0200521 using std::sqrt;
Benoit Jacobe2775862010-04-28 18:51:38 -0400522 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100523 typedef typename NumTraits<Scalar>::Real RealScalar;
524
Hauke Heibelf1679c72010-06-20 17:37:56 +0200525 Index rows = m.rows();
526 Index cols = m.cols();
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100527
Benoit Jacobe2775862010-04-28 18:51:38 -0400528 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200529 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200530 m3(rows, cols),
531 m4 = m1;
Benoit Steiner83149622015-12-10 13:13:45 -0800532
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800533 m4 = (m4.abs()==Scalar(0)).select(Scalar(1),m4);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100534
Gael Guennebaud9b1ad5e2012-03-09 12:08:06 +0100535 Scalar s1 = internal::random<Scalar>();
536
Deanna Hood41b717d2015-03-18 03:11:03 +1000537 // these tests are mostly to check possible compilation issues with free-functions.
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100538 VERIFY_IS_APPROX(m1.sin(), sin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100539 VERIFY_IS_APPROX(m1.cos(), cos(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000540 VERIFY_IS_APPROX(m1.tan(), tan(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100541 VERIFY_IS_APPROX(m1.asin(), asin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100542 VERIFY_IS_APPROX(m1.acos(), acos(m1));
Jitse Niesende150b12014-06-19 15:12:33 +0100543 VERIFY_IS_APPROX(m1.atan(), atan(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000544 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
545 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
546 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Rasmus Munk Larsen1e1848f2022-09-28 20:46:49 +0000547 VERIFY_IS_APPROX(m1.atan2(m2), atan2(m1,m2));
548
Rasmus Munk Larsen28ba1b22019-01-11 17:45:37 -0800549 VERIFY_IS_APPROX(m1.tanh().atanh(), atanh(tanh(m1)));
550 VERIFY_IS_APPROX(m1.sinh().asinh(), asinh(sinh(m1)));
551 VERIFY_IS_APPROX(m1.cosh().acosh(), acosh(cosh(m1)));
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700552 VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
Gael Guennebaud2f7e2612016-07-08 11:13:55 +0200553
Deanna Hooda5e49972015-03-11 08:56:42 +1000554 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000555 VERIFY_IS_APPROX(m1.round(), round(m1));
Ilya Tokar19876ce2019-12-16 16:00:35 -0500556 VERIFY_IS_APPROX(m1.rint(), rint(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000557 VERIFY_IS_APPROX(m1.floor(), floor(m1));
558 VERIFY_IS_APPROX(m1.ceil(), ceil(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200559 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
560 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
561 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800562 VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
Deanna Hood41b717d2015-03-18 03:11:03 +1000563 VERIFY_IS_APPROX(m1.abs(), abs(m1));
564 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
565 VERIFY_IS_APPROX(m1.square(), square(m1));
566 VERIFY_IS_APPROX(m1.cube(), cube(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100567 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100568 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Rasmus Munk Larsen6964ae82020-07-07 01:54:04 +0000569 VERIFY((m1.sqrt().sign().isNaN() == (Eigen::isnan)(sign(sqrt(m1)))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000570
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800571 // avoid inf and NaNs so verification doesn't fail
572 m3 = m4.abs();
573 VERIFY_IS_APPROX(m3.sqrt(), sqrt(abs(m3)));
574 VERIFY_IS_APPROX(m3.rsqrt(), Scalar(1)/sqrt(abs(m3)));
575 VERIFY_IS_APPROX(rsqrt(m3), Scalar(1)/sqrt(abs(m3)));
Deanna Hood41b717d2015-03-18 03:11:03 +1000576 VERIFY_IS_APPROX(m3.log(), log(m3));
Gael Guennebaud89099b02016-06-01 17:00:08 +0200577 VERIFY_IS_APPROX(m3.log1p(), log1p(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000578 VERIFY_IS_APPROX(m3.log10(), log10(m3));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000579 VERIFY_IS_APPROX(m3.log2(), log2(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000580
581
582 VERIFY((!(m1>m2) == (m1<=m2)).all());
583
584 VERIFY_IS_APPROX(sin(m1.asin()), m1);
585 VERIFY_IS_APPROX(cos(m1.acos()), m1);
586 VERIFY_IS_APPROX(tan(m1.atan()), m1);
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800587 VERIFY_IS_APPROX(sinh(m1), Scalar(0.5)*(exp(m1)-exp(-m1)));
588 VERIFY_IS_APPROX(cosh(m1), Scalar(0.5)*(exp(m1)+exp(-m1)));
589 VERIFY_IS_APPROX(tanh(m1), (Scalar(0.5)*(exp(m1)-exp(-m1)))/(Scalar(0.5)*(exp(m1)+exp(-m1))));
590 VERIFY_IS_APPROX(logistic(m1), (Scalar(1)/(Scalar(1)+exp(-m1))));
591 VERIFY_IS_APPROX(arg(m1), ((m1<Scalar(0)).template cast<Scalar>())*Scalar(std::acos(Scalar(-1))));
Deanna Hood41b717d2015-03-18 03:11:03 +1000592 VERIFY((round(m1) <= ceil(m1) && round(m1) >= floor(m1)).all());
Ilya Tokar19876ce2019-12-16 16:00:35 -0500593 VERIFY((rint(m1) <= ceil(m1) && rint(m1) >= floor(m1)).all());
594 VERIFY(((ceil(m1) - round(m1)) <= Scalar(0.5) || (round(m1) - floor(m1)) <= Scalar(0.5)).all());
595 VERIFY(((ceil(m1) - round(m1)) <= Scalar(1.0) && (round(m1) - floor(m1)) <= Scalar(1.0)).all());
596 VERIFY(((ceil(m1) - rint(m1)) <= Scalar(0.5) || (rint(m1) - floor(m1)) <= Scalar(0.5)).all());
597 VERIFY(((ceil(m1) - rint(m1)) <= Scalar(1.0) && (rint(m1) - floor(m1)) <= Scalar(1.0)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800598 VERIFY((Eigen::isnan)((m1*Scalar(0))/Scalar(0)).all());
599 VERIFY((Eigen::isinf)(m4/Scalar(0)).all());
600 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*Scalar(0)/Scalar(0))) && (!(Eigen::isfinite)(m4/Scalar(0)))).all());
601 VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
Deanna Hood41b717d2015-03-18 03:11:03 +1000602 VERIFY((abs(m1) == m1 || abs(m1) == -m1).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800603 VERIFY_IS_APPROX(m3, sqrt(abs2(m3)));
Joel Holdsworthd5c66572020-03-19 17:45:20 +0000604 VERIFY_IS_APPROX(m1.absolute_difference(m2), (m1 > m2).select(m1 - m2, m2 - m1));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100605 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
606 VERIFY_IS_APPROX( m1*m1.sign(),m1.abs());
607 VERIFY_IS_APPROX(m1.sign() * m1.abs(), m1);
Rasmus Munk Larsen1e1848f2022-09-28 20:46:49 +0000608
609 ArrayType tmp = m1.atan2(m2);
610 for (Index i = 0; i < tmp.size(); ++i) {
611 Scalar actual = tmp.array()(i);
612 Scalar expected = atan2(m1.array()(i), m2.array()(i));
613 VERIFY_IS_APPROX(actual, expected);
614 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400615
Gael Guennebaud62670c82013-06-10 23:40:56 +0200616 VERIFY_IS_APPROX(numext::abs2(numext::real(m1)) + numext::abs2(numext::imag(m1)), numext::abs2(m1));
Gael Guennebaud87427d22019-10-08 09:15:17 +0200617 VERIFY_IS_APPROX(numext::abs2(Eigen::real(m1)) + numext::abs2(Eigen::imag(m1)), numext::abs2(m1));
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400618 if(!NumTraits<Scalar>::IsComplex)
Gael Guennebaud62670c82013-06-10 23:40:56 +0200619 VERIFY_IS_APPROX(numext::real(m1), m1);
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200620
Jitse Niesen6fecb6f2014-02-24 14:10:17 +0000621 // shift argument of logarithm so that it is not zero
622 Scalar smallNumber = NumTraits<Scalar>::dummy_precision();
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800623 VERIFY_IS_APPROX((m3 + smallNumber).log() , log(abs(m3) + smallNumber));
624 VERIFY_IS_APPROX((m3 + smallNumber + Scalar(1)).log() , log1p(abs(m3) + smallNumber));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200625
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100626 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
627 VERIFY_IS_APPROX(m1.exp(), exp(m1));
628 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
Gael Guennebaud575ac542010-06-19 23:17:07 +0200629
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800630 VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800631 VERIFY_IS_APPROX((m3 + smallNumber).exp() - Scalar(1), expm1(abs(m3) + smallNumber));
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800632
Gael Guennebaud575ac542010-06-19 23:17:07 +0200633 VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100634 VERIFY_IS_APPROX(pow(m3,RealScalar(0.5)), m3.sqrt());
Benoit Steinere25e3a02015-12-03 18:16:35 -0800635
636 VERIFY_IS_APPROX(m3.pow(RealScalar(-0.5)), m3.rsqrt());
637 VERIFY_IS_APPROX(pow(m3,RealScalar(-0.5)), m3.rsqrt());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800638
639 // Avoid inf and NaN.
640 m3 = (m1.square()<NumTraits<Scalar>::epsilon()).select(Scalar(1),m3);
641 VERIFY_IS_APPROX(m3.pow(RealScalar(-2)), m3.square().inverse());
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +0000642
643 // Test pow and atan2 on special IEEE values.
644 binary_ops_test<Scalar>();
645 pow_scalar_exponent_test<Scalar>();
Benoit Steinere25e3a02015-12-03 18:16:35 -0800646
Antonio Sanchez4c42d5e2021-01-23 11:54:00 -0800647 VERIFY_IS_APPROX(log10(m3), log(m3)/numext::log(Scalar(10)));
648 VERIFY_IS_APPROX(log2(m3), log(m3)/numext::log(Scalar(2)));
Hauke Heibel81c13362012-03-07 08:58:42 +0100649
650 // scalar by array division
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100651 const RealScalar tiny = sqrt(std::numeric_limits<RealScalar>::epsilon());
Hauke Heibeldd9365e2012-03-09 14:04:13 +0100652 s1 += Scalar(tiny);
653 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
Rasmus Munk Larsen96dc37a2022-01-07 01:10:17 +0000654 VERIFY_IS_CWISE_APPROX(s1/m1, s1 * m1.inverse());
Eugene Brevdof7362772015-12-24 21:15:38 -0800655
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200656 // check inplace transpose
657 m3 = m1;
658 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000659 VERIFY_IS_APPROX(m3, m1.transpose());
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200660 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000661 VERIFY_IS_APPROX(m3, m1);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100662}
663
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000664
Jitse Niesen837db082011-05-09 10:17:41 +0100665template<typename ArrayType> void array_complex(const ArrayType& m)
666{
Deanna Hood41b717d2015-03-18 03:11:03 +1000667 typedef typename ArrayType::Scalar Scalar;
668 typedef typename NumTraits<Scalar>::Real RealScalar;
Jitse Niesen837db082011-05-09 10:17:41 +0100669
670 Index rows = m.rows();
671 Index cols = m.cols();
672
673 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200674 m2(rows, cols),
675 m4 = m1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800676
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200677 m4.real() = (m4.real().abs()==RealScalar(0)).select(RealScalar(1),m4.real());
678 m4.imag() = (m4.imag().abs()==RealScalar(0)).select(RealScalar(1),m4.imag());
Jitse Niesen837db082011-05-09 10:17:41 +0100679
Deanna Hood41b717d2015-03-18 03:11:03 +1000680 Array<RealScalar, -1, -1> m3(rows, cols);
681
Jitse Niesen837db082011-05-09 10:17:41 +0100682 for (Index i = 0; i < m.rows(); ++i)
683 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100684 m2(i,j) = sqrt(m1(i,j));
Jitse Niesen837db082011-05-09 10:17:41 +0100685
Deanna Hood41b717d2015-03-18 03:11:03 +1000686 // these tests are mostly to check possible compilation issues with free-functions.
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000687 VERIFY_IS_APPROX(m1.sin(), sin(m1));
688 VERIFY_IS_APPROX(m1.cos(), cos(m1));
689 VERIFY_IS_APPROX(m1.tan(), tan(m1));
690 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
691 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
692 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700693 VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000694 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200695 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
696 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
697 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800698 VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
Deanna Hood41b717d2015-03-18 03:11:03 +1000699 VERIFY_IS_APPROX(m1.log(), log(m1));
700 VERIFY_IS_APPROX(m1.log10(), log10(m1));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000701 VERIFY_IS_APPROX(m1.log2(), log2(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000702 VERIFY_IS_APPROX(m1.abs(), abs(m1));
703 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
704 VERIFY_IS_APPROX(m1.sqrt(), sqrt(m1));
705 VERIFY_IS_APPROX(m1.square(), square(m1));
706 VERIFY_IS_APPROX(m1.cube(), cube(m1));
707 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100708 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000709
Deanna Hood41b717d2015-03-18 03:11:03 +1000710 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
711 VERIFY_IS_APPROX(m1.exp(), exp(m1));
712 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
713
Srinivas Vasudevan88062b72019-12-12 01:56:54 +0000714 VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
715 VERIFY_IS_APPROX(expm1(m1), exp(m1) - 1.);
716 // Check for larger magnitude complex numbers that expm1 matches exp - 1.
717 VERIFY_IS_APPROX(expm1(10. * m1), exp(10. * m1) - 1.);
718
Deanna Hood41b717d2015-03-18 03:11:03 +1000719 VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
720 VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
721 VERIFY_IS_APPROX(tanh(m1), (0.5*(exp(m1)-exp(-m1)))/(0.5*(exp(m1)+exp(-m1))));
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700722 VERIFY_IS_APPROX(logistic(m1), (1.0/(1.0 + exp(-m1))));
Deanna Hood41b717d2015-03-18 03:11:03 +1000723
724 for (Index i = 0; i < m.rows(); ++i)
725 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebaud87427d22019-10-08 09:15:17 +0200726 m3(i,j) = std::atan2(m1(i,j).imag(), m1(i,j).real());
Deanna Hood41b717d2015-03-18 03:11:03 +1000727 VERIFY_IS_APPROX(arg(m1), m3);
728
729 std::complex<RealScalar> zero(0.0,0.0);
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200730 VERIFY((Eigen::isnan)(m1*zero/zero).all());
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100731#if EIGEN_COMP_MSVC
732 // msvc complex division is not robust
733 VERIFY((Eigen::isinf)(m4/RealScalar(0)).all());
734#else
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200735#if EIGEN_COMP_CLANG
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100736 // clang's complex division is notoriously broken too
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200737 if((numext::isinf)(m4(0,0)/RealScalar(0))) {
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200738#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100739 VERIFY((Eigen::isinf)(m4/zero).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200740#if EIGEN_COMP_CLANG
741 }
742 else
743 {
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200744 VERIFY((Eigen::isinf)(m4.real()/zero.real()).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200745 }
746#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100747#endif // MSVC
748
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200749 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*zero/zero)) && (!(Eigen::isfinite)(m1/zero))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000750
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800751 VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
Deanna Hood41b717d2015-03-18 03:11:03 +1000752 VERIFY_IS_APPROX(conj(m1.conjugate()), m1);
Gael Guennebaud87427d22019-10-08 09:15:17 +0200753 VERIFY_IS_APPROX(abs(m1), sqrt(square(m1.real())+square(m1.imag())));
Deanna Hood41b717d2015-03-18 03:11:03 +1000754 VERIFY_IS_APPROX(abs(m1), sqrt(abs2(m1)));
755 VERIFY_IS_APPROX(log10(m1), log(m1)/log(10));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000756 VERIFY_IS_APPROX(log2(m1), log(m1)/log(2));
Deanna Hood41b717d2015-03-18 03:11:03 +1000757
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100758 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
759 VERIFY_IS_APPROX( m1.sign() * m1.abs(), m1);
760
Deanna Hood41b717d2015-03-18 03:11:03 +1000761 // scalar by array division
Eugene Brevdof7362772015-12-24 21:15:38 -0800762 Scalar s1 = internal::random<Scalar>();
Benoit Steiner2d74ef92016-05-17 07:20:11 -0700763 const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon());
Deanna Hood41b717d2015-03-18 03:11:03 +1000764 s1 += Scalar(tiny);
765 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
766 VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
767
768 // check inplace transpose
769 m2 = m1;
770 m2.transposeInPlace();
771 VERIFY_IS_APPROX(m2, m1.transpose());
772 m2.transposeInPlace();
773 VERIFY_IS_APPROX(m2, m1);
Rasmus Munk Larsenb47c7772020-04-27 18:55:15 -0700774 // Check vectorized inplace transpose.
Rasmus Munk Larsen74ec8e62020-05-07 17:29:56 +0000775 ArrayType m5 = ArrayType::Random(131, 131);
Rasmus Munk Larsenb47c7772020-04-27 18:55:15 -0700776 ArrayType m6 = m5;
777 m6.transposeInPlace();
778 VERIFY_IS_APPROX(m6, m5.transpose());
Jitse Niesen837db082011-05-09 10:17:41 +0100779}
780
Abraham Bachrach039408c2012-01-11 11:00:30 -0500781template<typename ArrayType> void min_max(const ArrayType& m)
782{
Abraham Bachrach039408c2012-01-11 11:00:30 -0500783 typedef typename ArrayType::Scalar Scalar;
784
785 Index rows = m.rows();
786 Index cols = m.cols();
787
788 ArrayType m1 = ArrayType::Random(rows, cols);
789
790 // min/max with array
791 Scalar maxM1 = m1.maxCoeff();
792 Scalar minM1 = m1.minCoeff();
793
794 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1)));
795 VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1)));
796
797 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1)));
798 VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1)));
799
800 // min/max with scalar input
801 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1));
802 VERIFY_IS_APPROX(m1, (m1.min)( maxM1));
803
804 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1));
805 VERIFY_IS_APPROX(m1, (m1.max)( minM1));
806
Rasmus Munk Larsen841c8982021-02-24 17:49:20 -0800807
808 // min/max with various NaN propagation options.
809 if (m1.size() > 1 && !NumTraits<Scalar>::IsInteger) {
Antonio Sanchez8dfe1022021-03-16 20:12:46 -0700810 m1(0,0) = NumTraits<Scalar>::quiet_NaN();
Rasmus Munk Larsen841c8982021-02-24 17:49:20 -0800811 maxM1 = m1.template maxCoeff<PropagateNaN>();
812 minM1 = m1.template minCoeff<PropagateNaN>();
813 VERIFY((numext::isnan)(maxM1));
814 VERIFY((numext::isnan)(minM1));
815
816 maxM1 = m1.template maxCoeff<PropagateNumbers>();
817 minM1 = m1.template minCoeff<PropagateNumbers>();
818 VERIFY(!(numext::isnan)(maxM1));
819 VERIFY(!(numext::isnan)(minM1));
820 }
Abraham Bachrach039408c2012-01-11 11:00:30 -0500821}
822
Antonio Sanchezfc9d3522021-08-17 20:04:48 -0700823template<int N>
824struct shift_left {
825 template<typename Scalar>
826 Scalar operator()(const Scalar& v) const {
827 return v << N;
828 }
829};
830
831template<int N>
832struct arithmetic_shift_right {
833 template<typename Scalar>
834 Scalar operator()(const Scalar& v) const {
835 return v >> N;
836 }
837};
838
839template<typename ArrayType> void array_integer(const ArrayType& m)
840{
841 Index rows = m.rows();
842 Index cols = m.cols();
843
844 ArrayType m1 = ArrayType::Random(rows, cols),
845 m2(rows, cols);
846
847 m2 = m1.template shiftLeft<2>();
848 VERIFY( (m2 == m1.unaryExpr(shift_left<2>())).all() );
849 m2 = m1.template shiftLeft<9>();
850 VERIFY( (m2 == m1.unaryExpr(shift_left<9>())).all() );
851
852 m2 = m1.template shiftRight<2>();
853 VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<2>())).all() );
854 m2 = m1.template shiftRight<9>();
855 VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<9>())).all() );
856}
857
Gael Guennebaude0f6d352018-09-20 18:07:32 +0200858EIGEN_DECLARE_TEST(array_cwise)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000859{
860 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100861 CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100862 CALL_SUBTEST_2( array(Array22f()) );
863 CALL_SUBTEST_3( array(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200864 CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
865 CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
866 CALL_SUBTEST_6( array(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Gael Guennebaud543529d2019-01-22 15:30:50 +0100867 CALL_SUBTEST_6( array(Array<Index,Dynamic,Dynamic>(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Antonio Sanchezfc9d3522021-08-17 20:04:48 -0700868 CALL_SUBTEST_6( array_integer(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
869 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))) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000870 }
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100871 for(int i = 0; i < g_repeat; i++) {
872 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100873 CALL_SUBTEST_2( comparisons(Array22f()) );
874 CALL_SUBTEST_3( comparisons(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200875 CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
876 CALL_SUBTEST_6( comparisons(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100877 }
878 for(int i = 0; i < g_repeat; i++) {
Abraham Bachrach039408c2012-01-11 11:00:30 -0500879 CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) );
880 CALL_SUBTEST_2( min_max(Array22f()) );
881 CALL_SUBTEST_3( min_max(Array44d()) );
882 CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
883 CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
884 }
885 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100886 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100887 CALL_SUBTEST_2( array_real(Array22f()) );
888 CALL_SUBTEST_3( array_real(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200889 CALL_SUBTEST_5( array_real(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800890 CALL_SUBTEST_7( array_real(Array<Eigen::half, 32, 32>()) );
891 CALL_SUBTEST_8( array_real(Array<Eigen::bfloat16, 32, 32>()) );
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100892 }
Jitse Niesen837db082011-05-09 10:17:41 +0100893 for(int i = 0; i < g_repeat; i++) {
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200894 CALL_SUBTEST_4( array_complex(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
Jitse Niesen837db082011-05-09 10:17:41 +0100895 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400896
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000897 for(int i = 0; i < g_repeat; i++) {
898 CALL_SUBTEST_6( int_pow_test() );
899 CALL_SUBTEST_7( mixed_pow_test() );
900 }
901
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100902 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
903 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
904 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
Gael Guennebaud64fcfd32016-06-14 11:26:57 +0200905 typedef CwiseUnaryOp<internal::scalar_abs_op<double>, ArrayXd > Xpr;
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100906 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
907 ArrayBase<Xpr>
908 >::value));
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000909}