blob: 94c94517bbfe80fbe7ff65c6afc54201e4ec4f5e [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 }
Charles Schlosser82b152d2022-11-04 00:31:20 +0000222}
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000223
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
Charles Schlosser82b152d2022-11-04 00:31:20 +0000263namespace Eigen {
264namespace internal {
265template <typename Scalar>
266struct test_signbit_op {
267 Scalar constexpr operator()(const Scalar& a) const { return numext::signbit(a); }
268 template <typename Packet>
269 inline Packet packetOp(const Packet& a) const {
270 return psignbit(a);
271 }
272};
273template <typename Scalar>
274struct functor_traits<test_signbit_op<Scalar>> {
275 enum { Cost = 1, PacketAccess = true }; //todo: define HasSignbit flag
276};
277} // namespace internal
278} // namespace Eigen
279
280template <typename T, bool IsInteger = NumTraits<T>::IsInteger>
281struct ref_signbit_func_impl {
282 static bool run(const T& x) { return std::signbit(x); }
283};
284template <typename T>
285struct ref_signbit_func_impl<T, true> {
286 // MSVC (perhaps others) does not have a std::signbit overload for integers
287 static bool run(const T& x) { return x < T(0); }
288};
289template <typename T>
290bool ref_signbit_func(const T& x) {
291 return ref_signbit_func_impl<T>::run(x);
292}
293
294template <typename Scalar>
295void signbit_test() {
296 Scalar true_mask;
297 std::memset(static_cast<void*>(&true_mask), 0xff, sizeof(Scalar));
298 Scalar false_mask;
299 std::memset(static_cast<void*>(&false_mask), 0x00, sizeof(Scalar));
300
301 const size_t size = 100 * internal::packet_traits<Scalar>::size;
302 ArrayX<Scalar> x(size), y(size);
303 x.setRandom();
304 std::vector<Scalar> special_vals = special_values<Scalar>();
305 for (size_t i = 0; i < special_vals.size(); i++) {
306 x(2 * i + 0) = special_vals[i];
307 x(2 * i + 1) = -special_vals[i];
308 }
309 y = x.unaryExpr(internal::test_signbit_op<Scalar>());
310
311 bool all_pass = true;
312 for (size_t i = 0; i < size; i++) {
313 const Scalar ref_val = ref_signbit_func(x(i)) ? true_mask : false_mask;
314 bool not_same = internal::predux_any(internal::bitwise_helper<Scalar>::bitwise_xor(ref_val, y(i)));
315 if (not_same) std::cout << "signbit(" << x(i) << ") != " << y(i) << "\n";
316 all_pass = all_pass && !not_same;
317 }
318
319 VERIFY(all_pass);
320}
321void signbit_tests() {
322 signbit_test<float>();
323 signbit_test<double>();
324 signbit_test<Eigen::half>();
325 signbit_test<Eigen::bfloat16>();
326
327 signbit_test<uint8_t>();
328 signbit_test<uint16_t>();
329 signbit_test<uint32_t>();
330 signbit_test<uint64_t>();
331
332 signbit_test<int8_t>();
333 signbit_test<int16_t>();
334 signbit_test<int32_t>();
335 signbit_test<int64_t>();
336}
337
Benoit Jacobe2775862010-04-28 18:51:38 -0400338template<typename ArrayType> void array(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000339{
Benoit Jacobe2775862010-04-28 18:51:38 -0400340 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200341 typedef typename ArrayType::RealScalar RealScalar;
Benoit Jacobe2775862010-04-28 18:51:38 -0400342 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
343 typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000344
Hauke Heibelf1679c72010-06-20 17:37:56 +0200345 Index rows = m.rows();
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800346 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000347
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000348 ArrayType m1 = ArrayType::Random(rows, cols);
349 if (NumTraits<RealScalar>::IsInteger && NumTraits<RealScalar>::IsSigned
350 && !NumTraits<Scalar>::IsComplex) {
351 // Here we cap the size of the values in m1 such that pow(3)/cube()
352 // doesn't overflow and result in undefined behavior. Notice that because
353 // pow(int, int) promotes its inputs and output to double (according to
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000354 // the C++ standard), we have to make sure that the result fits in 53 bits
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000355 // for int64,
356 RealScalar max_val =
357 numext::mini(RealScalar(std::cbrt(NumTraits<RealScalar>::highest())),
358 RealScalar(std::cbrt(1LL << 53)))/2;
359 m1.array() = (m1.abs().array() <= max_val).select(m1, Scalar(max_val));
360 }
361 ArrayType m2 = ArrayType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000362 m3(rows, cols);
Christoph Hertzberg3be9f5c2015-04-16 13:25:20 +0200363 ArrayType m4 = m1; // copy constructor
364 VERIFY_IS_APPROX(m1, m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000365
Gael Guennebaud627595a2009-06-10 11:20:30 +0200366 ColVectorType cv1 = ColVectorType::Random(rows);
367 RowVectorType rv1 = RowVectorType::Random(cols);
368
Benoit Jacob47160402010-10-25 10:15:22 -0400369 Scalar s1 = internal::random<Scalar>(),
Hauke Heibel8cb3e362012-03-02 16:27:27 +0100370 s2 = internal::random<Scalar>();
371
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000372 // scalar addition
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100373 VERIFY_IS_APPROX(m1 + s1, s1 + m1);
Benoit Jacobe2775862010-04-28 18:51:38 -0400374 VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100375 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
Benoit Jacobe2775862010-04-28 18:51:38 -0400376 VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
377 VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
378 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000379 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100380 m3 += s2;
381 VERIFY_IS_APPROX(m3, m1 + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000382 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100383 m3 -= s1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800384 VERIFY_IS_APPROX(m3, m1 - s1);
385
Hauke Heibelf578dc72010-12-16 17:34:13 +0100386 // scalar operators via Maps
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000387 m3 = m1; m4 = m1;
388 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
389 VERIFY_IS_APPROX(m4, m3 - m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800390
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000391 m3 = m1; m4 = m1;
392 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
393 VERIFY_IS_APPROX(m4, m3 + m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800394
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000395 m3 = m1; m4 = m1;
396 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
397 VERIFY_IS_APPROX(m4, m3 * m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800398
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000399 m3 = m1; m4 = m1;
Hauke Heibelf578dc72010-12-16 17:34:13 +0100400 m2 = ArrayType::Random(rows,cols);
401 m2 = (m2==0).select(1,m2);
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000402 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
403 VERIFY_IS_APPROX(m4, m3 / m2);
Gael Guennebaud05ad0832008-07-19 13:03:23 +0000404
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000405 // reductions
Gael Guennebaud42af5872013-02-23 22:58:14 +0100406 VERIFY_IS_APPROX(m1.abs().colwise().sum().sum(), m1.abs().sum());
407 VERIFY_IS_APPROX(m1.abs().rowwise().sum().sum(), m1.abs().sum());
408 using std::abs;
409 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.colwise().sum().sum() - m1.sum()), m1.abs().sum());
410 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.rowwise().sum().sum() - m1.sum()), m1.abs().sum());
Gael Guennebaud28e139a2013-02-23 23:06:45 +0100411 if (!internal::isMuchSmallerThan(abs(m1.sum() - (m1+m2).sum()), m1.abs().sum(), test_precision<Scalar>()))
Hauke Heibel83e3c452010-12-16 18:53:02 +0100412 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud66e99ab2016-06-06 15:11:41 +0200413 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +0200414
415 // vector-wise ops
416 m3 = m1;
417 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
418 m3 = m1;
419 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
420 m3 = m1;
421 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
422 m3 = m1;
423 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800424
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200425 // Conversion from scalar
426 VERIFY_IS_APPROX((m3 = s1), ArrayType::Constant(rows,cols,s1));
427 VERIFY_IS_APPROX((m3 = 1), ArrayType::Constant(rows,cols,1));
428 VERIFY_IS_APPROX((m3.topLeftCorner(rows,cols) = 1), ArrayType::Constant(rows,cols,1));
429 typedef Array<Scalar,
430 ArrayType::RowsAtCompileTime==Dynamic?2:ArrayType::RowsAtCompileTime,
431 ArrayType::ColsAtCompileTime==Dynamic?2:ArrayType::ColsAtCompileTime,
432 ArrayType::Options> FixedArrayType;
Gael Guennebaud543529d2019-01-22 15:30:50 +0100433 {
434 FixedArrayType f1(s1);
435 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
436 FixedArrayType f2(numext::real(s1));
437 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
438 FixedArrayType f3((int)100*numext::real(s1));
439 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
440 f1.setRandom();
441 FixedArrayType f4(f1.data());
442 VERIFY_IS_APPROX(f4, f1);
443 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100444 {
445 FixedArrayType f1{s1};
446 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
447 FixedArrayType f2{numext::real(s1)};
448 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
449 FixedArrayType f3{(int)100*numext::real(s1)};
450 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
451 f1.setRandom();
452 FixedArrayType f4{f1.data()};
453 VERIFY_IS_APPROX(f4, f1);
454 }
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800455
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200456 // pow
457 VERIFY_IS_APPROX(m1.pow(2), m1.square());
458 VERIFY_IS_APPROX(pow(m1,2), m1.square());
459 VERIFY_IS_APPROX(m1.pow(3), m1.cube());
460 VERIFY_IS_APPROX(pow(m1,3), m1.cube());
461 VERIFY_IS_APPROX((-m1).pow(3), -m1.cube());
462 VERIFY_IS_APPROX(pow(2*m1,3), 8*m1.cube());
463 ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2));
464 VERIFY_IS_APPROX(Eigen::pow(m1,exponents), m1.square());
465 VERIFY_IS_APPROX(m1.pow(exponents), m1.square());
466 VERIFY_IS_APPROX(Eigen::pow(2*m1,exponents), 4*m1.square());
467 VERIFY_IS_APPROX((2*m1).pow(exponents), 4*m1.square());
468 VERIFY_IS_APPROX(Eigen::pow(m1,2*exponents), m1.square().square());
469 VERIFY_IS_APPROX(m1.pow(2*exponents), m1.square().square());
Gael Guennebaude61cee72016-07-04 11:49:03 +0200470 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 +0200471
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200472 // Check possible conflicts with 1D ctor
473 typedef Array<Scalar, Dynamic, 1> OneDArrayType;
Gael Guennebaud543529d2019-01-22 15:30:50 +0100474 {
475 OneDArrayType o1(rows);
476 VERIFY(o1.size()==rows);
477 OneDArrayType o2(static_cast<int>(rows));
478 VERIFY(o2.size()==rows);
479 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100480 {
481 OneDArrayType o1{rows};
482 VERIFY(o1.size()==rows);
483 OneDArrayType o4{int(rows)};
484 VERIFY(o4.size()==rows);
485 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100486 // Check possible conflicts with 2D ctor
487 typedef Array<Scalar, Dynamic, Dynamic> TwoDArrayType;
488 typedef Array<Scalar, 2, 1> ArrayType2;
489 {
490 TwoDArrayType o1(rows,cols);
491 VERIFY(o1.rows()==rows);
492 VERIFY(o1.cols()==cols);
493 TwoDArrayType o2(static_cast<int>(rows),static_cast<int>(cols));
494 VERIFY(o2.rows()==rows);
495 VERIFY(o2.cols()==cols);
496
497 ArrayType2 o3(rows,cols);
498 VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
499 ArrayType2 o4(static_cast<int>(rows),static_cast<int>(cols));
500 VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
501 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100502 {
503 TwoDArrayType o1{rows,cols};
504 VERIFY(o1.rows()==rows);
505 VERIFY(o1.cols()==cols);
506 TwoDArrayType o2{int(rows),int(cols)};
507 VERIFY(o2.rows()==rows);
508 VERIFY(o2.cols()==cols);
509
510 ArrayType2 o3{rows,cols};
511 VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
512 ArrayType2 o4{int(rows),int(cols)};
513 VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
514 }
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000515}
516
Benoit Jacobe2775862010-04-28 18:51:38 -0400517template<typename ArrayType> void comparisons(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000518{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100519 using std::abs;
Benoit Jacobe2775862010-04-28 18:51:38 -0400520 typedef typename ArrayType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +0000521 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000522
Hauke Heibelf1679c72010-06-20 17:37:56 +0200523 Index rows = m.rows();
524 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000525
Benoit Jacob47160402010-10-25 10:15:22 -0400526 Index r = internal::random<Index>(0, rows-1),
527 c = internal::random<Index>(0, cols-1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000528
Benoit Jacobe2775862010-04-28 18:51:38 -0400529 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200530 m2 = ArrayType::Random(rows, cols),
531 m3(rows, cols),
532 m4 = m1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800533
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200534 m4 = (m4.abs()==Scalar(0)).select(1,m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000535
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100536 VERIFY(((m1 + Scalar(1)) > m1).all());
537 VERIFY(((m1 - Scalar(1)) < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000538 if (rows*cols>1)
539 {
540 m3 = m1;
541 m3(r,c) += 1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100542 VERIFY(! (m1 < m3).all() );
543 VERIFY(! (m1 > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000544 }
Christoph Hertzberg494fa992015-05-07 17:28:40 +0200545 VERIFY(!(m1 > m2 && m1 < m2).any());
546 VERIFY((m1 <= m2 || m1 >= m2).all());
Gael Guennebaud627595a2009-06-10 11:20:30 +0200547
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200548 // comparisons array to scalar
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100549 VERIFY( (m1 != (m1(r,c)+1) ).any() );
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200550 VERIFY( (m1 > (m1(r,c)-1) ).any() );
551 VERIFY( (m1 < (m1(r,c)+1) ).any() );
552 VERIFY( (m1 == m1(r,c) ).any() );
553
554 // comparisons scalar to array
555 VERIFY( ( (m1(r,c)+1) != m1).any() );
556 VERIFY( ( (m1(r,c)-1) < m1).any() );
557 VERIFY( ( (m1(r,c)+1) > m1).any() );
558 VERIFY( ( m1(r,c) == m1).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200559
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000560 // test Select
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100561 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
562 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
Gael Guennebaud9f795582009-12-16 19:18:40 +0100563 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000564 for (int j=0; j<cols; ++j)
565 for (int i=0; i<rows; ++i)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100566 m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
Benoit Jacobe2775862010-04-28 18:51:38 -0400567 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
568 .select(ArrayType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000569 // shorter versions:
Benoit Jacobe2775862010-04-28 18:51:38 -0400570 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000571 .select(0,m1), m3);
Benoit Jacobe2775862010-04-28 18:51:38 -0400572 VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000573 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000574 // even shorter version:
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100575 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200576
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000577 // count
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100578 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
Benoit Jacobaaaade42010-05-30 16:00:58 -0400579
Gael Guennebaud35c11582011-05-31 22:17:34 +0200580 // and/or
581 VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0);
582 VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols);
583 RealScalar a = m1.abs().mean();
584 VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count());
585
Gael Guennebaud12e1ebb2018-07-12 17:16:40 +0200586 typedef Array<Index, Dynamic, 1> ArrayOfIndices;
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200587
Gael Guennebaud9f795582009-12-16 19:18:40 +0100588 // TODO allows colwise/rowwise for array
Benoit Jacobaaaade42010-05-30 16:00:58 -0400589 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
590 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
Benoit Jacobe8009992008-11-03 22:47:00 +0000591}
592
Benoit Jacobe2775862010-04-28 18:51:38 -0400593template<typename ArrayType> void array_real(const ArrayType& m)
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100594{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100595 using std::abs;
Gael Guennebaud9b9177f2013-07-05 13:35:34 +0200596 using std::sqrt;
Benoit Jacobe2775862010-04-28 18:51:38 -0400597 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100598 typedef typename NumTraits<Scalar>::Real RealScalar;
599
Hauke Heibelf1679c72010-06-20 17:37:56 +0200600 Index rows = m.rows();
601 Index cols = m.cols();
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100602
Benoit Jacobe2775862010-04-28 18:51:38 -0400603 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200604 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200605 m3(rows, cols),
606 m4 = m1;
Benoit Steiner83149622015-12-10 13:13:45 -0800607
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800608 m4 = (m4.abs()==Scalar(0)).select(Scalar(1),m4);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100609
Gael Guennebaud9b1ad5e2012-03-09 12:08:06 +0100610 Scalar s1 = internal::random<Scalar>();
611
Deanna Hood41b717d2015-03-18 03:11:03 +1000612 // these tests are mostly to check possible compilation issues with free-functions.
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100613 VERIFY_IS_APPROX(m1.sin(), sin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100614 VERIFY_IS_APPROX(m1.cos(), cos(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000615 VERIFY_IS_APPROX(m1.tan(), tan(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100616 VERIFY_IS_APPROX(m1.asin(), asin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100617 VERIFY_IS_APPROX(m1.acos(), acos(m1));
Jitse Niesende150b12014-06-19 15:12:33 +0100618 VERIFY_IS_APPROX(m1.atan(), atan(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000619 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
620 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
621 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Rasmus Munk Larsen1e1848f2022-09-28 20:46:49 +0000622 VERIFY_IS_APPROX(m1.atan2(m2), atan2(m1,m2));
623
Rasmus Munk Larsen28ba1b22019-01-11 17:45:37 -0800624 VERIFY_IS_APPROX(m1.tanh().atanh(), atanh(tanh(m1)));
625 VERIFY_IS_APPROX(m1.sinh().asinh(), asinh(sinh(m1)));
626 VERIFY_IS_APPROX(m1.cosh().acosh(), acosh(cosh(m1)));
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700627 VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
Gael Guennebaud2f7e2612016-07-08 11:13:55 +0200628
Deanna Hooda5e49972015-03-11 08:56:42 +1000629 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000630 VERIFY_IS_APPROX(m1.round(), round(m1));
Ilya Tokar19876ce2019-12-16 16:00:35 -0500631 VERIFY_IS_APPROX(m1.rint(), rint(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000632 VERIFY_IS_APPROX(m1.floor(), floor(m1));
633 VERIFY_IS_APPROX(m1.ceil(), ceil(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200634 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
635 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
636 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800637 VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
Deanna Hood41b717d2015-03-18 03:11:03 +1000638 VERIFY_IS_APPROX(m1.abs(), abs(m1));
639 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
640 VERIFY_IS_APPROX(m1.square(), square(m1));
641 VERIFY_IS_APPROX(m1.cube(), cube(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100642 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100643 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Rasmus Munk Larsen6964ae82020-07-07 01:54:04 +0000644 VERIFY((m1.sqrt().sign().isNaN() == (Eigen::isnan)(sign(sqrt(m1)))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000645
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800646 // avoid inf and NaNs so verification doesn't fail
647 m3 = m4.abs();
648 VERIFY_IS_APPROX(m3.sqrt(), sqrt(abs(m3)));
649 VERIFY_IS_APPROX(m3.rsqrt(), Scalar(1)/sqrt(abs(m3)));
650 VERIFY_IS_APPROX(rsqrt(m3), Scalar(1)/sqrt(abs(m3)));
Deanna Hood41b717d2015-03-18 03:11:03 +1000651 VERIFY_IS_APPROX(m3.log(), log(m3));
Gael Guennebaud89099b02016-06-01 17:00:08 +0200652 VERIFY_IS_APPROX(m3.log1p(), log1p(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000653 VERIFY_IS_APPROX(m3.log10(), log10(m3));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000654 VERIFY_IS_APPROX(m3.log2(), log2(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000655
656
657 VERIFY((!(m1>m2) == (m1<=m2)).all());
658
659 VERIFY_IS_APPROX(sin(m1.asin()), m1);
660 VERIFY_IS_APPROX(cos(m1.acos()), m1);
661 VERIFY_IS_APPROX(tan(m1.atan()), m1);
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800662 VERIFY_IS_APPROX(sinh(m1), Scalar(0.5)*(exp(m1)-exp(-m1)));
663 VERIFY_IS_APPROX(cosh(m1), Scalar(0.5)*(exp(m1)+exp(-m1)));
664 VERIFY_IS_APPROX(tanh(m1), (Scalar(0.5)*(exp(m1)-exp(-m1)))/(Scalar(0.5)*(exp(m1)+exp(-m1))));
665 VERIFY_IS_APPROX(logistic(m1), (Scalar(1)/(Scalar(1)+exp(-m1))));
666 VERIFY_IS_APPROX(arg(m1), ((m1<Scalar(0)).template cast<Scalar>())*Scalar(std::acos(Scalar(-1))));
Deanna Hood41b717d2015-03-18 03:11:03 +1000667 VERIFY((round(m1) <= ceil(m1) && round(m1) >= floor(m1)).all());
Ilya Tokar19876ce2019-12-16 16:00:35 -0500668 VERIFY((rint(m1) <= ceil(m1) && rint(m1) >= floor(m1)).all());
669 VERIFY(((ceil(m1) - round(m1)) <= Scalar(0.5) || (round(m1) - floor(m1)) <= Scalar(0.5)).all());
670 VERIFY(((ceil(m1) - round(m1)) <= Scalar(1.0) && (round(m1) - floor(m1)) <= Scalar(1.0)).all());
671 VERIFY(((ceil(m1) - rint(m1)) <= Scalar(0.5) || (rint(m1) - floor(m1)) <= Scalar(0.5)).all());
672 VERIFY(((ceil(m1) - rint(m1)) <= Scalar(1.0) && (rint(m1) - floor(m1)) <= Scalar(1.0)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800673 VERIFY((Eigen::isnan)((m1*Scalar(0))/Scalar(0)).all());
674 VERIFY((Eigen::isinf)(m4/Scalar(0)).all());
675 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*Scalar(0)/Scalar(0))) && (!(Eigen::isfinite)(m4/Scalar(0)))).all());
676 VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
Deanna Hood41b717d2015-03-18 03:11:03 +1000677 VERIFY((abs(m1) == m1 || abs(m1) == -m1).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800678 VERIFY_IS_APPROX(m3, sqrt(abs2(m3)));
Joel Holdsworthd5c66572020-03-19 17:45:20 +0000679 VERIFY_IS_APPROX(m1.absolute_difference(m2), (m1 > m2).select(m1 - m2, m2 - m1));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100680 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
681 VERIFY_IS_APPROX( m1*m1.sign(),m1.abs());
682 VERIFY_IS_APPROX(m1.sign() * m1.abs(), m1);
Rasmus Munk Larsen1e1848f2022-09-28 20:46:49 +0000683
684 ArrayType tmp = m1.atan2(m2);
685 for (Index i = 0; i < tmp.size(); ++i) {
686 Scalar actual = tmp.array()(i);
687 Scalar expected = atan2(m1.array()(i), m2.array()(i));
688 VERIFY_IS_APPROX(actual, expected);
689 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400690
Gael Guennebaud62670c82013-06-10 23:40:56 +0200691 VERIFY_IS_APPROX(numext::abs2(numext::real(m1)) + numext::abs2(numext::imag(m1)), numext::abs2(m1));
Gael Guennebaud87427d22019-10-08 09:15:17 +0200692 VERIFY_IS_APPROX(numext::abs2(Eigen::real(m1)) + numext::abs2(Eigen::imag(m1)), numext::abs2(m1));
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400693 if(!NumTraits<Scalar>::IsComplex)
Gael Guennebaud62670c82013-06-10 23:40:56 +0200694 VERIFY_IS_APPROX(numext::real(m1), m1);
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200695
Jitse Niesen6fecb6f2014-02-24 14:10:17 +0000696 // shift argument of logarithm so that it is not zero
697 Scalar smallNumber = NumTraits<Scalar>::dummy_precision();
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800698 VERIFY_IS_APPROX((m3 + smallNumber).log() , log(abs(m3) + smallNumber));
699 VERIFY_IS_APPROX((m3 + smallNumber + Scalar(1)).log() , log1p(abs(m3) + smallNumber));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200700
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100701 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
702 VERIFY_IS_APPROX(m1.exp(), exp(m1));
703 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
Gael Guennebaud575ac542010-06-19 23:17:07 +0200704
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800705 VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800706 VERIFY_IS_APPROX((m3 + smallNumber).exp() - Scalar(1), expm1(abs(m3) + smallNumber));
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800707
Gael Guennebaud575ac542010-06-19 23:17:07 +0200708 VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100709 VERIFY_IS_APPROX(pow(m3,RealScalar(0.5)), m3.sqrt());
Benoit Steinere25e3a02015-12-03 18:16:35 -0800710
711 VERIFY_IS_APPROX(m3.pow(RealScalar(-0.5)), m3.rsqrt());
712 VERIFY_IS_APPROX(pow(m3,RealScalar(-0.5)), m3.rsqrt());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800713
714 // Avoid inf and NaN.
715 m3 = (m1.square()<NumTraits<Scalar>::epsilon()).select(Scalar(1),m3);
716 VERIFY_IS_APPROX(m3.pow(RealScalar(-2)), m3.square().inverse());
Rasmus Munk Larsen14c847d2022-10-12 20:12:08 +0000717
718 // Test pow and atan2 on special IEEE values.
719 binary_ops_test<Scalar>();
720 pow_scalar_exponent_test<Scalar>();
Benoit Steinere25e3a02015-12-03 18:16:35 -0800721
Antonio Sanchez4c42d5e2021-01-23 11:54:00 -0800722 VERIFY_IS_APPROX(log10(m3), log(m3)/numext::log(Scalar(10)));
723 VERIFY_IS_APPROX(log2(m3), log(m3)/numext::log(Scalar(2)));
Hauke Heibel81c13362012-03-07 08:58:42 +0100724
725 // scalar by array division
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100726 const RealScalar tiny = sqrt(std::numeric_limits<RealScalar>::epsilon());
Hauke Heibeldd9365e2012-03-09 14:04:13 +0100727 s1 += Scalar(tiny);
728 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
Rasmus Munk Larsen96dc37a2022-01-07 01:10:17 +0000729 VERIFY_IS_CWISE_APPROX(s1/m1, s1 * m1.inverse());
Eugene Brevdof7362772015-12-24 21:15:38 -0800730
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200731 // check inplace transpose
732 m3 = m1;
733 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000734 VERIFY_IS_APPROX(m3, m1.transpose());
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200735 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000736 VERIFY_IS_APPROX(m3, m1);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100737}
738
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000739
Jitse Niesen837db082011-05-09 10:17:41 +0100740template<typename ArrayType> void array_complex(const ArrayType& m)
741{
Deanna Hood41b717d2015-03-18 03:11:03 +1000742 typedef typename ArrayType::Scalar Scalar;
743 typedef typename NumTraits<Scalar>::Real RealScalar;
Jitse Niesen837db082011-05-09 10:17:41 +0100744
745 Index rows = m.rows();
746 Index cols = m.cols();
747
748 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200749 m2(rows, cols),
750 m4 = m1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800751
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200752 m4.real() = (m4.real().abs()==RealScalar(0)).select(RealScalar(1),m4.real());
753 m4.imag() = (m4.imag().abs()==RealScalar(0)).select(RealScalar(1),m4.imag());
Jitse Niesen837db082011-05-09 10:17:41 +0100754
Deanna Hood41b717d2015-03-18 03:11:03 +1000755 Array<RealScalar, -1, -1> m3(rows, cols);
756
Jitse Niesen837db082011-05-09 10:17:41 +0100757 for (Index i = 0; i < m.rows(); ++i)
758 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100759 m2(i,j) = sqrt(m1(i,j));
Jitse Niesen837db082011-05-09 10:17:41 +0100760
Deanna Hood41b717d2015-03-18 03:11:03 +1000761 // these tests are mostly to check possible compilation issues with free-functions.
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000762 VERIFY_IS_APPROX(m1.sin(), sin(m1));
763 VERIFY_IS_APPROX(m1.cos(), cos(m1));
764 VERIFY_IS_APPROX(m1.tan(), tan(m1));
765 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
766 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
767 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700768 VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000769 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200770 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
771 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
772 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800773 VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
Deanna Hood41b717d2015-03-18 03:11:03 +1000774 VERIFY_IS_APPROX(m1.log(), log(m1));
775 VERIFY_IS_APPROX(m1.log10(), log10(m1));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000776 VERIFY_IS_APPROX(m1.log2(), log2(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000777 VERIFY_IS_APPROX(m1.abs(), abs(m1));
778 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
779 VERIFY_IS_APPROX(m1.sqrt(), sqrt(m1));
780 VERIFY_IS_APPROX(m1.square(), square(m1));
781 VERIFY_IS_APPROX(m1.cube(), cube(m1));
782 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100783 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000784
Deanna Hood41b717d2015-03-18 03:11:03 +1000785 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
786 VERIFY_IS_APPROX(m1.exp(), exp(m1));
787 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
788
Srinivas Vasudevan88062b72019-12-12 01:56:54 +0000789 VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
790 VERIFY_IS_APPROX(expm1(m1), exp(m1) - 1.);
791 // Check for larger magnitude complex numbers that expm1 matches exp - 1.
792 VERIFY_IS_APPROX(expm1(10. * m1), exp(10. * m1) - 1.);
793
Deanna Hood41b717d2015-03-18 03:11:03 +1000794 VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
795 VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
796 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 -0700797 VERIFY_IS_APPROX(logistic(m1), (1.0/(1.0 + exp(-m1))));
Deanna Hood41b717d2015-03-18 03:11:03 +1000798
799 for (Index i = 0; i < m.rows(); ++i)
800 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebaud87427d22019-10-08 09:15:17 +0200801 m3(i,j) = std::atan2(m1(i,j).imag(), m1(i,j).real());
Deanna Hood41b717d2015-03-18 03:11:03 +1000802 VERIFY_IS_APPROX(arg(m1), m3);
803
804 std::complex<RealScalar> zero(0.0,0.0);
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200805 VERIFY((Eigen::isnan)(m1*zero/zero).all());
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100806#if EIGEN_COMP_MSVC
807 // msvc complex division is not robust
808 VERIFY((Eigen::isinf)(m4/RealScalar(0)).all());
809#else
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200810#if EIGEN_COMP_CLANG
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100811 // clang's complex division is notoriously broken too
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200812 if((numext::isinf)(m4(0,0)/RealScalar(0))) {
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200813#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100814 VERIFY((Eigen::isinf)(m4/zero).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200815#if EIGEN_COMP_CLANG
816 }
817 else
818 {
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200819 VERIFY((Eigen::isinf)(m4.real()/zero.real()).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200820 }
821#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100822#endif // MSVC
823
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200824 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*zero/zero)) && (!(Eigen::isfinite)(m1/zero))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000825
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800826 VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
Deanna Hood41b717d2015-03-18 03:11:03 +1000827 VERIFY_IS_APPROX(conj(m1.conjugate()), m1);
Gael Guennebaud87427d22019-10-08 09:15:17 +0200828 VERIFY_IS_APPROX(abs(m1), sqrt(square(m1.real())+square(m1.imag())));
Deanna Hood41b717d2015-03-18 03:11:03 +1000829 VERIFY_IS_APPROX(abs(m1), sqrt(abs2(m1)));
830 VERIFY_IS_APPROX(log10(m1), log(m1)/log(10));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000831 VERIFY_IS_APPROX(log2(m1), log(m1)/log(2));
Deanna Hood41b717d2015-03-18 03:11:03 +1000832
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100833 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
834 VERIFY_IS_APPROX( m1.sign() * m1.abs(), m1);
835
Deanna Hood41b717d2015-03-18 03:11:03 +1000836 // scalar by array division
Eugene Brevdof7362772015-12-24 21:15:38 -0800837 Scalar s1 = internal::random<Scalar>();
Benoit Steiner2d74ef92016-05-17 07:20:11 -0700838 const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon());
Deanna Hood41b717d2015-03-18 03:11:03 +1000839 s1 += Scalar(tiny);
840 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
841 VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
842
843 // check inplace transpose
844 m2 = m1;
845 m2.transposeInPlace();
846 VERIFY_IS_APPROX(m2, m1.transpose());
847 m2.transposeInPlace();
848 VERIFY_IS_APPROX(m2, m1);
Rasmus Munk Larsenb47c7772020-04-27 18:55:15 -0700849 // Check vectorized inplace transpose.
Rasmus Munk Larsen74ec8e62020-05-07 17:29:56 +0000850 ArrayType m5 = ArrayType::Random(131, 131);
Rasmus Munk Larsenb47c7772020-04-27 18:55:15 -0700851 ArrayType m6 = m5;
852 m6.transposeInPlace();
853 VERIFY_IS_APPROX(m6, m5.transpose());
Jitse Niesen837db082011-05-09 10:17:41 +0100854}
855
Abraham Bachrach039408c2012-01-11 11:00:30 -0500856template<typename ArrayType> void min_max(const ArrayType& m)
857{
Abraham Bachrach039408c2012-01-11 11:00:30 -0500858 typedef typename ArrayType::Scalar Scalar;
859
860 Index rows = m.rows();
861 Index cols = m.cols();
862
863 ArrayType m1 = ArrayType::Random(rows, cols);
864
865 // min/max with array
866 Scalar maxM1 = m1.maxCoeff();
867 Scalar minM1 = m1.minCoeff();
868
869 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1)));
870 VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1)));
871
872 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1)));
873 VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1)));
874
875 // min/max with scalar input
876 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1));
877 VERIFY_IS_APPROX(m1, (m1.min)( maxM1));
878
879 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1));
880 VERIFY_IS_APPROX(m1, (m1.max)( minM1));
881
Rasmus Munk Larsen841c8982021-02-24 17:49:20 -0800882
883 // min/max with various NaN propagation options.
884 if (m1.size() > 1 && !NumTraits<Scalar>::IsInteger) {
Antonio Sanchez8dfe1022021-03-16 20:12:46 -0700885 m1(0,0) = NumTraits<Scalar>::quiet_NaN();
Rasmus Munk Larsen841c8982021-02-24 17:49:20 -0800886 maxM1 = m1.template maxCoeff<PropagateNaN>();
887 minM1 = m1.template minCoeff<PropagateNaN>();
888 VERIFY((numext::isnan)(maxM1));
889 VERIFY((numext::isnan)(minM1));
890
891 maxM1 = m1.template maxCoeff<PropagateNumbers>();
892 minM1 = m1.template minCoeff<PropagateNumbers>();
893 VERIFY(!(numext::isnan)(maxM1));
894 VERIFY(!(numext::isnan)(minM1));
895 }
Abraham Bachrach039408c2012-01-11 11:00:30 -0500896}
897
Antonio Sanchezfc9d3522021-08-17 20:04:48 -0700898template<int N>
899struct shift_left {
900 template<typename Scalar>
901 Scalar operator()(const Scalar& v) const {
902 return v << N;
903 }
904};
905
906template<int N>
907struct arithmetic_shift_right {
908 template<typename Scalar>
909 Scalar operator()(const Scalar& v) const {
910 return v >> N;
911 }
912};
913
914template<typename ArrayType> void array_integer(const ArrayType& m)
915{
916 Index rows = m.rows();
917 Index cols = m.cols();
918
919 ArrayType m1 = ArrayType::Random(rows, cols),
920 m2(rows, cols);
921
922 m2 = m1.template shiftLeft<2>();
923 VERIFY( (m2 == m1.unaryExpr(shift_left<2>())).all() );
924 m2 = m1.template shiftLeft<9>();
925 VERIFY( (m2 == m1.unaryExpr(shift_left<9>())).all() );
926
927 m2 = m1.template shiftRight<2>();
928 VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<2>())).all() );
929 m2 = m1.template shiftRight<9>();
930 VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<9>())).all() );
931}
932
Charles Schlosser82b152d2022-11-04 00:31:20 +0000933template <typename ArrayType>
934struct signed_shift_test_impl {
935 typedef typename ArrayType::Scalar Scalar;
936 static constexpr size_t Size = sizeof(Scalar);
937 static constexpr size_t MaxShift = (CHAR_BIT * Size) - 1;
938
939 template <size_t N = 0>
940 static inline std::enable_if_t<(N > MaxShift), void> run(const ArrayType& ) {}
941 template <size_t N = 0>
942 static inline std::enable_if_t<(N <= MaxShift), void> run(const ArrayType& m) {
943 const Index rows = m.rows();
944 const Index cols = m.cols();
945
946 ArrayType m1 = ArrayType::Random(rows, cols), m2(rows, cols);
947
948 m2 = m1.unaryExpr([](const Scalar& x) { return x >> N; });
949 VERIFY((m2 == m1.unaryExpr(internal::scalar_shift_right_op<Scalar, N>())).all());
950
951 m2 = m1.unaryExpr([](const Scalar& x) { return x << N; });
952 VERIFY((m2 == m1.unaryExpr( internal::scalar_shift_left_op<Scalar, N>())).all());
953
954 run<N + 1>(m);
955 }
956};
957template <typename ArrayType>
958void signed_shift_test(const ArrayType& m) {
959 signed_shift_test_impl<ArrayType>::run(m);
960}
961
Gael Guennebaude0f6d352018-09-20 18:07:32 +0200962EIGEN_DECLARE_TEST(array_cwise)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000963{
964 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100965 CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100966 CALL_SUBTEST_2( array(Array22f()) );
967 CALL_SUBTEST_3( array(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200968 CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
969 CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
970 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 +0100971 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 -0700972 CALL_SUBTEST_6( array_integer(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
973 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))) );
Charles Schlosser82b152d2022-11-04 00:31:20 +0000974 CALL_SUBTEST_7( signed_shift_test(ArrayXXi(internal::random<int>(1, EIGEN_TEST_MAX_SIZE), internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
975 CALL_SUBTEST_7( signed_shift_test(Array<Index, Dynamic, Dynamic>(internal::random<int>(1, EIGEN_TEST_MAX_SIZE), internal::random<int>(1, EIGEN_TEST_MAX_SIZE))));
976
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000977 }
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100978 for(int i = 0; i < g_repeat; i++) {
979 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100980 CALL_SUBTEST_2( comparisons(Array22f()) );
981 CALL_SUBTEST_3( comparisons(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200982 CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
983 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 +0100984 }
985 for(int i = 0; i < g_repeat; i++) {
Abraham Bachrach039408c2012-01-11 11:00:30 -0500986 CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) );
987 CALL_SUBTEST_2( min_max(Array22f()) );
988 CALL_SUBTEST_3( min_max(Array44d()) );
989 CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
990 CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
991 }
992 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100993 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100994 CALL_SUBTEST_2( array_real(Array22f()) );
995 CALL_SUBTEST_3( array_real(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200996 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 -0800997 CALL_SUBTEST_7( array_real(Array<Eigen::half, 32, 32>()) );
998 CALL_SUBTEST_8( array_real(Array<Eigen::bfloat16, 32, 32>()) );
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100999 }
Jitse Niesen837db082011-05-09 10:17:41 +01001000 for(int i = 0; i < g_repeat; i++) {
Gael Guennebauda8f66fe2011-07-12 14:41:00 +02001001 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 +01001002 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -04001003
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +00001004 for(int i = 0; i < g_repeat; i++) {
1005 CALL_SUBTEST_6( int_pow_test() );
1006 CALL_SUBTEST_7( mixed_pow_test() );
Charles Schlosser82b152d2022-11-04 00:31:20 +00001007 CALL_SUBTEST_8( signbit_tests() );
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +00001008 }
1009
Hauke Heibel2d0dfe52010-12-16 17:36:10 +01001010 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
1011 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
1012 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
Gael Guennebaud64fcfd32016-06-14 11:26:57 +02001013 typedef CwiseUnaryOp<internal::scalar_abs_op<double>, ArrayXd > Xpr;
Hauke Heibel2d0dfe52010-12-16 17:36:10 +01001014 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
1015 ArrayBase<Xpr>
1016 >::value));
Gael Guennebaud8cef5412008-06-21 17:28:07 +00001017}