blob: 5a4c282782d1bace85deed3c65e73443f2f08ff0 [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
Gael Guennebaud8cef5412008-06-21 17:28:07 +000010#include "main.h"
Gael Guennebaud8cef5412008-06-21 17:28:07 +000011
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000012
13// Test the corner cases of pow(x, y) for real types.
14template<typename Scalar>
15void pow_test() {
16 const Scalar zero = Scalar(0);
Antonio Sanchez8dfe1022021-03-16 20:12:46 -070017 const Scalar eps = Eigen::NumTraits<Scalar>::epsilon();
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000018 const Scalar one = Scalar(1);
Rasmus Munk Larsen6e3b7952021-02-05 16:58:49 -080019 const Scalar two = Scalar(2);
20 const Scalar three = Scalar(3);
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000021 const Scalar sqrt_half = Scalar(std::sqrt(0.5));
22 const Scalar sqrt2 = Scalar(std::sqrt(2));
Antonio Sanchez8dfe1022021-03-16 20:12:46 -070023 const Scalar inf = Eigen::NumTraits<Scalar>::infinity();
24 const Scalar nan = Eigen::NumTraits<Scalar>::quiet_NaN();
Rasmus Munk Larsenbe0574e2021-02-17 02:50:32 +000025 const Scalar denorm_min = std::numeric_limits<Scalar>::denorm_min();
Rasmus Munk Larsen6e3b7952021-02-05 16:58:49 -080026 const Scalar min = (std::numeric_limits<Scalar>::min)();
27 const Scalar max = (std::numeric_limits<Scalar>::max)();
Antonio Sanchez8dfe1022021-03-16 20:12:46 -070028 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 +000029
Rasmus Munk Larsen6e3b7952021-02-05 16:58:49 -080030 const static Scalar abs_vals[] = {zero,
Rasmus Munk Larsenbe0574e2021-02-17 02:50:32 +000031 denorm_min,
32 min,
33 eps,
Rasmus Munk Larsen6e3b7952021-02-05 16:58:49 -080034 sqrt_half,
35 one,
36 sqrt2,
37 two,
38 three,
Rasmus Munk Larsenbe0574e2021-02-17 02:50:32 +000039 max_exp,
Rasmus Munk Larsen6e3b7952021-02-05 16:58:49 -080040 max,
41 inf,
42 nan};
Rasmus Munk Larsenbe0574e2021-02-17 02:50:32 +000043 const int abs_cases = 13;
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000044 const int num_cases = 2*abs_cases * 2*abs_cases;
45 // Repeat the same value to make sure we hit the vectorized path.
46 const int num_repeats = 32;
47 Array<Scalar, Dynamic, Dynamic> x(num_repeats, num_cases);
48 Array<Scalar, Dynamic, Dynamic> y(num_repeats, num_cases);
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000049 int count = 0;
50 for (int i = 0; i < abs_cases; ++i) {
51 const Scalar abs_x = abs_vals[i];
52 for (int sign_x = 0; sign_x < 2; ++sign_x) {
53 Scalar x_case = sign_x == 0 ? -abs_x : abs_x;
54 for (int j = 0; j < abs_cases; ++j) {
55 const Scalar abs_y = abs_vals[j];
56 for (int sign_y = 0; sign_y < 2; ++sign_y) {
57 Scalar y_case = sign_y == 0 ? -abs_y : abs_y;
58 for (int repeat = 0; repeat < num_repeats; ++repeat) {
59 x(repeat, count) = x_case;
60 y(repeat, count) = y_case;
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +000061 }
62 ++count;
63 }
64 }
65 }
66 }
67
68 Array<Scalar, Dynamic, Dynamic> actual = x.pow(y);
69 const Scalar tol = test_precision<Scalar>();
70 bool all_pass = true;
71 for (int i = 0; i < 1; ++i) {
72 for (int j = 0; j < num_cases; ++j) {
Rasmus Munk Larsenbe0574e2021-02-17 02:50:32 +000073 Scalar e = static_cast<Scalar>(std::pow(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 Larsencdd8fdc2021-01-18 13:25:16 +000078 std::cout << "pow(" << x(i,j) << "," << y(i,j) << ") = " << a << " != " << e << std::endl;
79 }
80 }
81 }
Charles Schlosser76a669f2022-08-16 21:32:36 +000082
83 typedef typename internal::make_integer<Scalar>::type Int_t;
84
85 // ensure both vectorized and non-vectorized paths taken
86 Index test_size = 2 * internal::packet_traits<Scalar>::size + 1;
87
88 Array<Scalar, Dynamic, 1> eigenPow(test_size);
89 for (int i = 0; i < num_cases; ++i) {
90 Array<Scalar, Dynamic, 1> bases = x.col(i);
91 for (Scalar abs_exponent : abs_vals){
92 for (Scalar exponent : {-abs_exponent, abs_exponent}){
93 // test floating point exponent code path
94 eigenPow.setZero();
95 eigenPow = bases.pow(exponent);
96 for (int j = 0; j < num_repeats; j++){
97 Scalar e = static_cast<Scalar>(std::pow(bases(j), exponent));
98 Scalar a = eigenPow(j);
99 bool success = (a == e) || ((numext::isfinite)(e) && internal::isApprox(a, e, tol)) || ((numext::isnan)(a) && (numext::isnan)(e));
100 all_pass &= success;
101 if (!success) {
102 std::cout << "pow(" << x(i, j) << "," << y(i, j) << ") = " << a << " != " << e << std::endl;
103 }
104 }
105 // test integer exponent code path
106 bool exponent_is_integer = (numext::isfinite)(exponent) && (numext::round(exponent) == exponent) && (numext::abs(exponent) < static_cast<Scalar>(NumTraits<Int_t>::highest()));
107 if (exponent_is_integer)
108 {
109 Int_t exponent_as_int = static_cast<Int_t>(exponent);
110 eigenPow.setZero();
111 eigenPow = bases.pow(exponent_as_int);
112 for (int j = 0; j < num_repeats; j++){
113 Scalar e = static_cast<Scalar>(std::pow(bases(j), exponent));
114 Scalar a = eigenPow(j);
115 bool success = (a == e) || ((numext::isfinite)(e) && internal::isApprox(a, e, tol)) || ((numext::isnan)(a) && (numext::isnan)(e));
116 all_pass &= success;
117 if (!success) {
118 std::cout << "pow(" << x(i, j) << "," << y(i, j) << ") = " << a << " != " << e << std::endl;
119 }
120 }
121 }
122 }
123 }
124 }
125
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +0000126 VERIFY(all_pass);
127}
128
Charles Schlossere5af9f82022-08-29 19:23:54 +0000129template <typename Scalar, typename ScalarExponent>
130Scalar calc_overflow_threshold(const ScalarExponent exponent) {
131 EIGEN_USING_STD(exp2);
Antonio Sánchez69f50e32022-09-06 19:53:29 +0000132 EIGEN_USING_STD(log2);
Charles Schlossere5af9f82022-08-29 19:23:54 +0000133 EIGEN_STATIC_ASSERT((NumTraits<Scalar>::digits() < 2 * NumTraits<double>::digits()), BASE_TYPE_IS_TOO_BIG);
134
135 if (exponent < 2)
136 return NumTraits<Scalar>::highest();
137 else {
Antonio Sánchez69f50e32022-09-06 19:53:29 +0000138 // base^e <= highest ==> base <= 2^(log2(highest)/e)
Antonio Sánchez3c37dd22022-09-08 17:40:45 +0000139 // For floating-point types, consider the bound for integer values that can be reproduced exactly = 2 ^ digits
140 double highest_bits = numext::mini(static_cast<double>(NumTraits<Scalar>::digits()),
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000141 static_cast<double>(log2(NumTraits<Scalar>::highest())));
Antonio Sánchez69f50e32022-09-06 19:53:29 +0000142 return static_cast<Scalar>(
Antonio Sánchez3c37dd22022-09-08 17:40:45 +0000143 numext::floor(exp2(highest_bits / static_cast<double>(exponent))));
Charles Schlossere5af9f82022-08-29 19:23:54 +0000144 }
145}
146
147template <typename Base, typename Exponent>
148void test_exponent(Exponent exponent) {
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000149 const Base max_abs_bases = static_cast<Base>(10000);
150 // avoid integer overflow in Base type
151 Base threshold = calc_overflow_threshold<Base, Exponent>(numext::abs(exponent));
152 // avoid numbers that can't be verified with std::pow
153 double double_threshold = calc_overflow_threshold<double, Exponent>(numext::abs(exponent));
154 // use the lesser of these two thresholds
155 Base testing_threshold =
156 static_cast<double>(threshold) < double_threshold ? threshold : static_cast<Base>(double_threshold);
157 // test both vectorized and non-vectorized code paths
158 const Index array_size = 2 * internal::packet_traits<Base>::size + 1;
159
160 Base max_base = numext::mini(testing_threshold, max_abs_bases);
161 Base min_base = NumTraits<Base>::IsSigned ? -max_base : Base(0);
162
163 ArrayX<Base> x(array_size), y(array_size);
164 bool all_pass = true;
165 for (Base base = min_base; base <= max_base; base++) {
166 if (exponent < 0 && base == 0) continue;
167 x.setConstant(base);
168 y = x.pow(exponent);
Charles Schlossere5af9f82022-08-29 19:23:54 +0000169 EIGEN_USING_STD(pow);
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000170 Base e = pow(base, static_cast<Base>(exponent));
171 for (Base a : y) {
172 bool pass = (a == e);
173 if (!NumTraits<Base>::IsInteger) {
174 pass = pass || (((numext::isfinite)(e) && internal::isApprox(a, e)) ||
175 ((numext::isnan)(a) && (numext::isnan)(e)));
176 }
177 all_pass &= pass;
178 if (!pass) {
179 std::cout << "pow(" << base << "," << exponent << ") = " << a << " != " << e << std::endl;
180 }
Charles Schlossere5af9f82022-08-29 19:23:54 +0000181 }
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000182 }
183 VERIFY(all_pass);
Charles Schlossere5af9f82022-08-29 19:23:54 +0000184}
Charles Schlossere5af9f82022-08-29 19:23:54 +0000185
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000186template <typename Base, typename Exponent>
187void unary_pow_test() {
188 Exponent max_exponent = static_cast<Exponent>(NumTraits<Base>::digits());
189 Exponent min_exponent = static_cast<Exponent>(NumTraits<Exponent>::IsSigned ? -max_exponent : 0);
190
191 for (Exponent exponent = min_exponent; exponent < max_exponent; ++exponent) {
192 test_exponent<Base, Exponent>(exponent);
193 }
194};
195
196void mixed_pow_test() {
197 // The following cases will test promoting a smaller exponent type
198 // to a wider base type.
199 unary_pow_test<double, int>();
200 unary_pow_test<double, float>();
201 unary_pow_test<float, half>();
202 unary_pow_test<double, half>();
203 unary_pow_test<float, bfloat16>();
204 unary_pow_test<double, bfloat16>();
205
206 // Although in the following cases the exponent cannot be represented exactly
207 // in the base type, we do not perform a conversion, but implement
208 // the operation using repeated squaring.
209 unary_pow_test<float, int>();
210 unary_pow_test<double, long long>();
211
212 // The following cases will test promoting a wider exponent type
213 // to a narrower base type. This should compile but generate a
214 // deprecation warning:
215 unary_pow_test<float, double>();
216}
217
218void int_pow_test() {
219 unary_pow_test<int, int>();
220 unary_pow_test<unsigned int, unsigned int>();
221 unary_pow_test<long long, long long>();
222 unary_pow_test<unsigned long long, unsigned long long>();
223
224 // Although in the following cases the exponent cannot be represented exactly
225 // in the base type, we do not perform a conversion, but implement the
226 // operation using repeated squaring.
227 unary_pow_test<long long, int>();
228 unary_pow_test<int, unsigned int>();
229 unary_pow_test<unsigned int, int>();
230 unary_pow_test<long long, unsigned long long>();
231 unary_pow_test<unsigned long long, long long>();
232 unary_pow_test<long long, int>();
Charles Schlossere5af9f82022-08-29 19:23:54 +0000233}
234
Benoit Jacobe2775862010-04-28 18:51:38 -0400235template<typename ArrayType> void array(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000236{
Benoit Jacobe2775862010-04-28 18:51:38 -0400237 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200238 typedef typename ArrayType::RealScalar RealScalar;
Benoit Jacobe2775862010-04-28 18:51:38 -0400239 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
240 typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000241
Hauke Heibelf1679c72010-06-20 17:37:56 +0200242 Index rows = m.rows();
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800243 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000244
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000245 ArrayType m1 = ArrayType::Random(rows, cols);
246 if (NumTraits<RealScalar>::IsInteger && NumTraits<RealScalar>::IsSigned
247 && !NumTraits<Scalar>::IsComplex) {
248 // Here we cap the size of the values in m1 such that pow(3)/cube()
249 // doesn't overflow and result in undefined behavior. Notice that because
250 // pow(int, int) promotes its inputs and output to double (according to
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000251 // the C++ standard), we have to make sure that the result fits in 53 bits
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000252 // for int64,
253 RealScalar max_val =
254 numext::mini(RealScalar(std::cbrt(NumTraits<RealScalar>::highest())),
255 RealScalar(std::cbrt(1LL << 53)))/2;
256 m1.array() = (m1.abs().array() <= max_val).select(m1, Scalar(max_val));
257 }
258 ArrayType m2 = ArrayType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000259 m3(rows, cols);
Christoph Hertzberg3be9f5c2015-04-16 13:25:20 +0200260 ArrayType m4 = m1; // copy constructor
261 VERIFY_IS_APPROX(m1, m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000262
Gael Guennebaud627595a2009-06-10 11:20:30 +0200263 ColVectorType cv1 = ColVectorType::Random(rows);
264 RowVectorType rv1 = RowVectorType::Random(cols);
265
Benoit Jacob47160402010-10-25 10:15:22 -0400266 Scalar s1 = internal::random<Scalar>(),
Hauke Heibel8cb3e362012-03-02 16:27:27 +0100267 s2 = internal::random<Scalar>();
268
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000269 // scalar addition
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100270 VERIFY_IS_APPROX(m1 + s1, s1 + m1);
Benoit Jacobe2775862010-04-28 18:51:38 -0400271 VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100272 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
Benoit Jacobe2775862010-04-28 18:51:38 -0400273 VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
274 VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
275 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000276 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100277 m3 += s2;
278 VERIFY_IS_APPROX(m3, m1 + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000279 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100280 m3 -= s1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800281 VERIFY_IS_APPROX(m3, m1 - s1);
282
Hauke Heibelf578dc72010-12-16 17:34:13 +0100283 // scalar operators via Maps
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000284 m3 = m1; m4 = m1;
285 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
286 VERIFY_IS_APPROX(m4, m3 - m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800287
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000288 m3 = m1; m4 = m1;
289 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
290 VERIFY_IS_APPROX(m4, m3 + m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800291
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000292 m3 = m1; m4 = m1;
293 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
294 VERIFY_IS_APPROX(m4, m3 * m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800295
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000296 m3 = m1; m4 = m1;
Hauke Heibelf578dc72010-12-16 17:34:13 +0100297 m2 = ArrayType::Random(rows,cols);
298 m2 = (m2==0).select(1,m2);
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000299 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
300 VERIFY_IS_APPROX(m4, m3 / m2);
Gael Guennebaud05ad0832008-07-19 13:03:23 +0000301
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000302 // reductions
Gael Guennebaud42af5872013-02-23 22:58:14 +0100303 VERIFY_IS_APPROX(m1.abs().colwise().sum().sum(), m1.abs().sum());
304 VERIFY_IS_APPROX(m1.abs().rowwise().sum().sum(), m1.abs().sum());
305 using std::abs;
306 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.colwise().sum().sum() - m1.sum()), m1.abs().sum());
307 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.rowwise().sum().sum() - m1.sum()), m1.abs().sum());
Gael Guennebaud28e139a2013-02-23 23:06:45 +0100308 if (!internal::isMuchSmallerThan(abs(m1.sum() - (m1+m2).sum()), m1.abs().sum(), test_precision<Scalar>()))
Hauke Heibel83e3c452010-12-16 18:53:02 +0100309 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud66e99ab2016-06-06 15:11:41 +0200310 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +0200311
312 // vector-wise ops
313 m3 = m1;
314 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
315 m3 = m1;
316 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
317 m3 = m1;
318 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
319 m3 = m1;
320 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800321
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200322 // Conversion from scalar
323 VERIFY_IS_APPROX((m3 = s1), ArrayType::Constant(rows,cols,s1));
324 VERIFY_IS_APPROX((m3 = 1), ArrayType::Constant(rows,cols,1));
325 VERIFY_IS_APPROX((m3.topLeftCorner(rows,cols) = 1), ArrayType::Constant(rows,cols,1));
326 typedef Array<Scalar,
327 ArrayType::RowsAtCompileTime==Dynamic?2:ArrayType::RowsAtCompileTime,
328 ArrayType::ColsAtCompileTime==Dynamic?2:ArrayType::ColsAtCompileTime,
329 ArrayType::Options> FixedArrayType;
Gael Guennebaud543529d2019-01-22 15:30:50 +0100330 {
331 FixedArrayType f1(s1);
332 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
333 FixedArrayType f2(numext::real(s1));
334 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
335 FixedArrayType f3((int)100*numext::real(s1));
336 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
337 f1.setRandom();
338 FixedArrayType f4(f1.data());
339 VERIFY_IS_APPROX(f4, f1);
340 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100341 {
342 FixedArrayType f1{s1};
343 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
344 FixedArrayType f2{numext::real(s1)};
345 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
346 FixedArrayType f3{(int)100*numext::real(s1)};
347 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
348 f1.setRandom();
349 FixedArrayType f4{f1.data()};
350 VERIFY_IS_APPROX(f4, f1);
351 }
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800352
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200353 // pow
354 VERIFY_IS_APPROX(m1.pow(2), m1.square());
355 VERIFY_IS_APPROX(pow(m1,2), m1.square());
356 VERIFY_IS_APPROX(m1.pow(3), m1.cube());
357 VERIFY_IS_APPROX(pow(m1,3), m1.cube());
358 VERIFY_IS_APPROX((-m1).pow(3), -m1.cube());
359 VERIFY_IS_APPROX(pow(2*m1,3), 8*m1.cube());
360 ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2));
361 VERIFY_IS_APPROX(Eigen::pow(m1,exponents), m1.square());
362 VERIFY_IS_APPROX(m1.pow(exponents), m1.square());
363 VERIFY_IS_APPROX(Eigen::pow(2*m1,exponents), 4*m1.square());
364 VERIFY_IS_APPROX((2*m1).pow(exponents), 4*m1.square());
365 VERIFY_IS_APPROX(Eigen::pow(m1,2*exponents), m1.square().square());
366 VERIFY_IS_APPROX(m1.pow(2*exponents), m1.square().square());
Gael Guennebaude61cee72016-07-04 11:49:03 +0200367 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 +0200368
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200369 // Check possible conflicts with 1D ctor
370 typedef Array<Scalar, Dynamic, 1> OneDArrayType;
Gael Guennebaud543529d2019-01-22 15:30:50 +0100371 {
372 OneDArrayType o1(rows);
373 VERIFY(o1.size()==rows);
374 OneDArrayType o2(static_cast<int>(rows));
375 VERIFY(o2.size()==rows);
376 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100377 {
378 OneDArrayType o1{rows};
379 VERIFY(o1.size()==rows);
380 OneDArrayType o4{int(rows)};
381 VERIFY(o4.size()==rows);
382 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100383 // Check possible conflicts with 2D ctor
384 typedef Array<Scalar, Dynamic, Dynamic> TwoDArrayType;
385 typedef Array<Scalar, 2, 1> ArrayType2;
386 {
387 TwoDArrayType o1(rows,cols);
388 VERIFY(o1.rows()==rows);
389 VERIFY(o1.cols()==cols);
390 TwoDArrayType o2(static_cast<int>(rows),static_cast<int>(cols));
391 VERIFY(o2.rows()==rows);
392 VERIFY(o2.cols()==cols);
393
394 ArrayType2 o3(rows,cols);
395 VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
396 ArrayType2 o4(static_cast<int>(rows),static_cast<int>(cols));
397 VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
398 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100399 {
400 TwoDArrayType o1{rows,cols};
401 VERIFY(o1.rows()==rows);
402 VERIFY(o1.cols()==cols);
403 TwoDArrayType o2{int(rows),int(cols)};
404 VERIFY(o2.rows()==rows);
405 VERIFY(o2.cols()==cols);
406
407 ArrayType2 o3{rows,cols};
408 VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
409 ArrayType2 o4{int(rows),int(cols)};
410 VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
411 }
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000412}
413
Benoit Jacobe2775862010-04-28 18:51:38 -0400414template<typename ArrayType> void comparisons(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000415{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100416 using std::abs;
Benoit Jacobe2775862010-04-28 18:51:38 -0400417 typedef typename ArrayType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +0000418 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000419
Hauke Heibelf1679c72010-06-20 17:37:56 +0200420 Index rows = m.rows();
421 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000422
Benoit Jacob47160402010-10-25 10:15:22 -0400423 Index r = internal::random<Index>(0, rows-1),
424 c = internal::random<Index>(0, cols-1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000425
Benoit Jacobe2775862010-04-28 18:51:38 -0400426 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200427 m2 = ArrayType::Random(rows, cols),
428 m3(rows, cols),
429 m4 = m1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800430
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200431 m4 = (m4.abs()==Scalar(0)).select(1,m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000432
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100433 VERIFY(((m1 + Scalar(1)) > m1).all());
434 VERIFY(((m1 - Scalar(1)) < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000435 if (rows*cols>1)
436 {
437 m3 = m1;
438 m3(r,c) += 1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100439 VERIFY(! (m1 < m3).all() );
440 VERIFY(! (m1 > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000441 }
Christoph Hertzberg494fa992015-05-07 17:28:40 +0200442 VERIFY(!(m1 > m2 && m1 < m2).any());
443 VERIFY((m1 <= m2 || m1 >= m2).all());
Gael Guennebaud627595a2009-06-10 11:20:30 +0200444
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200445 // comparisons array to scalar
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100446 VERIFY( (m1 != (m1(r,c)+1) ).any() );
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200447 VERIFY( (m1 > (m1(r,c)-1) ).any() );
448 VERIFY( (m1 < (m1(r,c)+1) ).any() );
449 VERIFY( (m1 == m1(r,c) ).any() );
450
451 // comparisons scalar to array
452 VERIFY( ( (m1(r,c)+1) != m1).any() );
453 VERIFY( ( (m1(r,c)-1) < m1).any() );
454 VERIFY( ( (m1(r,c)+1) > m1).any() );
455 VERIFY( ( m1(r,c) == m1).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200456
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000457 // test Select
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100458 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
459 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
Gael Guennebaud9f795582009-12-16 19:18:40 +0100460 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000461 for (int j=0; j<cols; ++j)
462 for (int i=0; i<rows; ++i)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100463 m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
Benoit Jacobe2775862010-04-28 18:51:38 -0400464 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
465 .select(ArrayType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000466 // shorter versions:
Benoit Jacobe2775862010-04-28 18:51:38 -0400467 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000468 .select(0,m1), m3);
Benoit Jacobe2775862010-04-28 18:51:38 -0400469 VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000470 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000471 // even shorter version:
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100472 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200473
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000474 // count
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100475 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
Benoit Jacobaaaade42010-05-30 16:00:58 -0400476
Gael Guennebaud35c11582011-05-31 22:17:34 +0200477 // and/or
478 VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0);
479 VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols);
480 RealScalar a = m1.abs().mean();
481 VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count());
482
Gael Guennebaud12e1ebb2018-07-12 17:16:40 +0200483 typedef Array<Index, Dynamic, 1> ArrayOfIndices;
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200484
Gael Guennebaud9f795582009-12-16 19:18:40 +0100485 // TODO allows colwise/rowwise for array
Benoit Jacobaaaade42010-05-30 16:00:58 -0400486 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
487 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
Benoit Jacobe8009992008-11-03 22:47:00 +0000488}
489
Benoit Jacobe2775862010-04-28 18:51:38 -0400490template<typename ArrayType> void array_real(const ArrayType& m)
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100491{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100492 using std::abs;
Gael Guennebaud9b9177f2013-07-05 13:35:34 +0200493 using std::sqrt;
Benoit Jacobe2775862010-04-28 18:51:38 -0400494 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100495 typedef typename NumTraits<Scalar>::Real RealScalar;
496
Hauke Heibelf1679c72010-06-20 17:37:56 +0200497 Index rows = m.rows();
498 Index cols = m.cols();
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100499
Benoit Jacobe2775862010-04-28 18:51:38 -0400500 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200501 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200502 m3(rows, cols),
503 m4 = m1;
Benoit Steiner83149622015-12-10 13:13:45 -0800504
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800505 m4 = (m4.abs()==Scalar(0)).select(Scalar(1),m4);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100506
Gael Guennebaud9b1ad5e2012-03-09 12:08:06 +0100507 Scalar s1 = internal::random<Scalar>();
508
Deanna Hood41b717d2015-03-18 03:11:03 +1000509 // these tests are mostly to check possible compilation issues with free-functions.
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100510 VERIFY_IS_APPROX(m1.sin(), sin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100511 VERIFY_IS_APPROX(m1.cos(), cos(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000512 VERIFY_IS_APPROX(m1.tan(), tan(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100513 VERIFY_IS_APPROX(m1.asin(), asin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100514 VERIFY_IS_APPROX(m1.acos(), acos(m1));
Jitse Niesende150b12014-06-19 15:12:33 +0100515 VERIFY_IS_APPROX(m1.atan(), atan(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000516 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
517 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
518 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Rasmus Munk Larsen28ba1b22019-01-11 17:45:37 -0800519#if EIGEN_HAS_CXX11_MATH
520 VERIFY_IS_APPROX(m1.tanh().atanh(), atanh(tanh(m1)));
521 VERIFY_IS_APPROX(m1.sinh().asinh(), asinh(sinh(m1)));
522 VERIFY_IS_APPROX(m1.cosh().acosh(), acosh(cosh(m1)));
523#endif
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700524 VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
Gael Guennebaud2f7e2612016-07-08 11:13:55 +0200525
Deanna Hooda5e49972015-03-11 08:56:42 +1000526 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000527 VERIFY_IS_APPROX(m1.round(), round(m1));
Ilya Tokar19876ce2019-12-16 16:00:35 -0500528 VERIFY_IS_APPROX(m1.rint(), rint(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000529 VERIFY_IS_APPROX(m1.floor(), floor(m1));
530 VERIFY_IS_APPROX(m1.ceil(), ceil(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200531 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
532 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
533 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800534 VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
Deanna Hood41b717d2015-03-18 03:11:03 +1000535 VERIFY_IS_APPROX(m1.abs(), abs(m1));
536 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
537 VERIFY_IS_APPROX(m1.square(), square(m1));
538 VERIFY_IS_APPROX(m1.cube(), cube(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100539 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100540 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Rasmus Munk Larsen6964ae82020-07-07 01:54:04 +0000541 VERIFY((m1.sqrt().sign().isNaN() == (Eigen::isnan)(sign(sqrt(m1)))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000542
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800543 // avoid inf and NaNs so verification doesn't fail
544 m3 = m4.abs();
545 VERIFY_IS_APPROX(m3.sqrt(), sqrt(abs(m3)));
546 VERIFY_IS_APPROX(m3.rsqrt(), Scalar(1)/sqrt(abs(m3)));
547 VERIFY_IS_APPROX(rsqrt(m3), Scalar(1)/sqrt(abs(m3)));
Deanna Hood41b717d2015-03-18 03:11:03 +1000548 VERIFY_IS_APPROX(m3.log(), log(m3));
Gael Guennebaud89099b02016-06-01 17:00:08 +0200549 VERIFY_IS_APPROX(m3.log1p(), log1p(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000550 VERIFY_IS_APPROX(m3.log10(), log10(m3));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000551 VERIFY_IS_APPROX(m3.log2(), log2(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000552
553
554 VERIFY((!(m1>m2) == (m1<=m2)).all());
555
556 VERIFY_IS_APPROX(sin(m1.asin()), m1);
557 VERIFY_IS_APPROX(cos(m1.acos()), m1);
558 VERIFY_IS_APPROX(tan(m1.atan()), m1);
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800559 VERIFY_IS_APPROX(sinh(m1), Scalar(0.5)*(exp(m1)-exp(-m1)));
560 VERIFY_IS_APPROX(cosh(m1), Scalar(0.5)*(exp(m1)+exp(-m1)));
561 VERIFY_IS_APPROX(tanh(m1), (Scalar(0.5)*(exp(m1)-exp(-m1)))/(Scalar(0.5)*(exp(m1)+exp(-m1))));
562 VERIFY_IS_APPROX(logistic(m1), (Scalar(1)/(Scalar(1)+exp(-m1))));
563 VERIFY_IS_APPROX(arg(m1), ((m1<Scalar(0)).template cast<Scalar>())*Scalar(std::acos(Scalar(-1))));
Deanna Hood41b717d2015-03-18 03:11:03 +1000564 VERIFY((round(m1) <= ceil(m1) && round(m1) >= floor(m1)).all());
Ilya Tokar19876ce2019-12-16 16:00:35 -0500565 VERIFY((rint(m1) <= ceil(m1) && rint(m1) >= floor(m1)).all());
566 VERIFY(((ceil(m1) - round(m1)) <= Scalar(0.5) || (round(m1) - floor(m1)) <= Scalar(0.5)).all());
567 VERIFY(((ceil(m1) - round(m1)) <= Scalar(1.0) && (round(m1) - floor(m1)) <= Scalar(1.0)).all());
568 VERIFY(((ceil(m1) - rint(m1)) <= Scalar(0.5) || (rint(m1) - floor(m1)) <= Scalar(0.5)).all());
569 VERIFY(((ceil(m1) - rint(m1)) <= Scalar(1.0) && (rint(m1) - floor(m1)) <= Scalar(1.0)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800570 VERIFY((Eigen::isnan)((m1*Scalar(0))/Scalar(0)).all());
571 VERIFY((Eigen::isinf)(m4/Scalar(0)).all());
572 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*Scalar(0)/Scalar(0))) && (!(Eigen::isfinite)(m4/Scalar(0)))).all());
573 VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
Deanna Hood41b717d2015-03-18 03:11:03 +1000574 VERIFY((abs(m1) == m1 || abs(m1) == -m1).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800575 VERIFY_IS_APPROX(m3, sqrt(abs2(m3)));
Joel Holdsworthd5c66572020-03-19 17:45:20 +0000576 VERIFY_IS_APPROX(m1.absolute_difference(m2), (m1 > m2).select(m1 - m2, m2 - m1));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100577 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
578 VERIFY_IS_APPROX( m1*m1.sign(),m1.abs());
579 VERIFY_IS_APPROX(m1.sign() * m1.abs(), m1);
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400580
Gael Guennebaud62670c82013-06-10 23:40:56 +0200581 VERIFY_IS_APPROX(numext::abs2(numext::real(m1)) + numext::abs2(numext::imag(m1)), numext::abs2(m1));
Gael Guennebaud87427d22019-10-08 09:15:17 +0200582 VERIFY_IS_APPROX(numext::abs2(Eigen::real(m1)) + numext::abs2(Eigen::imag(m1)), numext::abs2(m1));
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400583 if(!NumTraits<Scalar>::IsComplex)
Gael Guennebaud62670c82013-06-10 23:40:56 +0200584 VERIFY_IS_APPROX(numext::real(m1), m1);
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200585
Jitse Niesen6fecb6f2014-02-24 14:10:17 +0000586 // shift argument of logarithm so that it is not zero
587 Scalar smallNumber = NumTraits<Scalar>::dummy_precision();
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800588 VERIFY_IS_APPROX((m3 + smallNumber).log() , log(abs(m3) + smallNumber));
589 VERIFY_IS_APPROX((m3 + smallNumber + Scalar(1)).log() , log1p(abs(m3) + smallNumber));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200590
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100591 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
592 VERIFY_IS_APPROX(m1.exp(), exp(m1));
593 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
Gael Guennebaud575ac542010-06-19 23:17:07 +0200594
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800595 VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800596 VERIFY_IS_APPROX((m3 + smallNumber).exp() - Scalar(1), expm1(abs(m3) + smallNumber));
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800597
Gael Guennebaud575ac542010-06-19 23:17:07 +0200598 VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100599 VERIFY_IS_APPROX(pow(m3,RealScalar(0.5)), m3.sqrt());
Benoit Steinere25e3a02015-12-03 18:16:35 -0800600
601 VERIFY_IS_APPROX(m3.pow(RealScalar(-0.5)), m3.rsqrt());
602 VERIFY_IS_APPROX(pow(m3,RealScalar(-0.5)), m3.rsqrt());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800603
604 // Avoid inf and NaN.
605 m3 = (m1.square()<NumTraits<Scalar>::epsilon()).select(Scalar(1),m3);
606 VERIFY_IS_APPROX(m3.pow(RealScalar(-2)), m3.square().inverse());
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +0000607 pow_test<Scalar>();
Benoit Steinere25e3a02015-12-03 18:16:35 -0800608
Antonio Sanchez4c42d5e2021-01-23 11:54:00 -0800609 VERIFY_IS_APPROX(log10(m3), log(m3)/numext::log(Scalar(10)));
610 VERIFY_IS_APPROX(log2(m3), log(m3)/numext::log(Scalar(2)));
Hauke Heibel81c13362012-03-07 08:58:42 +0100611
612 // scalar by array division
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100613 const RealScalar tiny = sqrt(std::numeric_limits<RealScalar>::epsilon());
Hauke Heibeldd9365e2012-03-09 14:04:13 +0100614 s1 += Scalar(tiny);
615 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
Rasmus Munk Larsen96dc37a2022-01-07 01:10:17 +0000616 VERIFY_IS_CWISE_APPROX(s1/m1, s1 * m1.inverse());
Eugene Brevdof7362772015-12-24 21:15:38 -0800617
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200618 // check inplace transpose
619 m3 = m1;
620 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000621 VERIFY_IS_APPROX(m3, m1.transpose());
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200622 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000623 VERIFY_IS_APPROX(m3, m1);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100624}
625
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000626
Jitse Niesen837db082011-05-09 10:17:41 +0100627template<typename ArrayType> void array_complex(const ArrayType& m)
628{
Deanna Hood41b717d2015-03-18 03:11:03 +1000629 typedef typename ArrayType::Scalar Scalar;
630 typedef typename NumTraits<Scalar>::Real RealScalar;
Jitse Niesen837db082011-05-09 10:17:41 +0100631
632 Index rows = m.rows();
633 Index cols = m.cols();
634
635 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200636 m2(rows, cols),
637 m4 = m1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800638
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200639 m4.real() = (m4.real().abs()==RealScalar(0)).select(RealScalar(1),m4.real());
640 m4.imag() = (m4.imag().abs()==RealScalar(0)).select(RealScalar(1),m4.imag());
Jitse Niesen837db082011-05-09 10:17:41 +0100641
Deanna Hood41b717d2015-03-18 03:11:03 +1000642 Array<RealScalar, -1, -1> m3(rows, cols);
643
Jitse Niesen837db082011-05-09 10:17:41 +0100644 for (Index i = 0; i < m.rows(); ++i)
645 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100646 m2(i,j) = sqrt(m1(i,j));
Jitse Niesen837db082011-05-09 10:17:41 +0100647
Deanna Hood41b717d2015-03-18 03:11:03 +1000648 // these tests are mostly to check possible compilation issues with free-functions.
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000649 VERIFY_IS_APPROX(m1.sin(), sin(m1));
650 VERIFY_IS_APPROX(m1.cos(), cos(m1));
651 VERIFY_IS_APPROX(m1.tan(), tan(m1));
652 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
653 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
654 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700655 VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000656 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200657 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
658 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
659 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800660 VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
Deanna Hood41b717d2015-03-18 03:11:03 +1000661 VERIFY_IS_APPROX(m1.log(), log(m1));
662 VERIFY_IS_APPROX(m1.log10(), log10(m1));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000663 VERIFY_IS_APPROX(m1.log2(), log2(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000664 VERIFY_IS_APPROX(m1.abs(), abs(m1));
665 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
666 VERIFY_IS_APPROX(m1.sqrt(), sqrt(m1));
667 VERIFY_IS_APPROX(m1.square(), square(m1));
668 VERIFY_IS_APPROX(m1.cube(), cube(m1));
669 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100670 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000671
672
673 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
674 VERIFY_IS_APPROX(m1.exp(), exp(m1));
675 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
676
Srinivas Vasudevan88062b72019-12-12 01:56:54 +0000677 VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
678 VERIFY_IS_APPROX(expm1(m1), exp(m1) - 1.);
679 // Check for larger magnitude complex numbers that expm1 matches exp - 1.
680 VERIFY_IS_APPROX(expm1(10. * m1), exp(10. * m1) - 1.);
681
Deanna Hood41b717d2015-03-18 03:11:03 +1000682 VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
683 VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
684 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 -0700685 VERIFY_IS_APPROX(logistic(m1), (1.0/(1.0 + exp(-m1))));
Deanna Hood41b717d2015-03-18 03:11:03 +1000686
687 for (Index i = 0; i < m.rows(); ++i)
688 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebaud87427d22019-10-08 09:15:17 +0200689 m3(i,j) = std::atan2(m1(i,j).imag(), m1(i,j).real());
Deanna Hood41b717d2015-03-18 03:11:03 +1000690 VERIFY_IS_APPROX(arg(m1), m3);
691
692 std::complex<RealScalar> zero(0.0,0.0);
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200693 VERIFY((Eigen::isnan)(m1*zero/zero).all());
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100694#if EIGEN_COMP_MSVC
695 // msvc complex division is not robust
696 VERIFY((Eigen::isinf)(m4/RealScalar(0)).all());
697#else
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200698#if EIGEN_COMP_CLANG
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100699 // clang's complex division is notoriously broken too
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200700 if((numext::isinf)(m4(0,0)/RealScalar(0))) {
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200701#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100702 VERIFY((Eigen::isinf)(m4/zero).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200703#if EIGEN_COMP_CLANG
704 }
705 else
706 {
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200707 VERIFY((Eigen::isinf)(m4.real()/zero.real()).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200708 }
709#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100710#endif // MSVC
711
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200712 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*zero/zero)) && (!(Eigen::isfinite)(m1/zero))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000713
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800714 VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
Deanna Hood41b717d2015-03-18 03:11:03 +1000715 VERIFY_IS_APPROX(conj(m1.conjugate()), m1);
Gael Guennebaud87427d22019-10-08 09:15:17 +0200716 VERIFY_IS_APPROX(abs(m1), sqrt(square(m1.real())+square(m1.imag())));
Deanna Hood41b717d2015-03-18 03:11:03 +1000717 VERIFY_IS_APPROX(abs(m1), sqrt(abs2(m1)));
718 VERIFY_IS_APPROX(log10(m1), log(m1)/log(10));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000719 VERIFY_IS_APPROX(log2(m1), log(m1)/log(2));
Deanna Hood41b717d2015-03-18 03:11:03 +1000720
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100721 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
722 VERIFY_IS_APPROX( m1.sign() * m1.abs(), m1);
723
Deanna Hood41b717d2015-03-18 03:11:03 +1000724 // scalar by array division
Eugene Brevdof7362772015-12-24 21:15:38 -0800725 Scalar s1 = internal::random<Scalar>();
Benoit Steiner2d74ef92016-05-17 07:20:11 -0700726 const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon());
Deanna Hood41b717d2015-03-18 03:11:03 +1000727 s1 += Scalar(tiny);
728 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
729 VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
730
731 // check inplace transpose
732 m2 = m1;
733 m2.transposeInPlace();
734 VERIFY_IS_APPROX(m2, m1.transpose());
735 m2.transposeInPlace();
736 VERIFY_IS_APPROX(m2, m1);
Rasmus Munk Larsenb47c7772020-04-27 18:55:15 -0700737 // Check vectorized inplace transpose.
Rasmus Munk Larsen74ec8e62020-05-07 17:29:56 +0000738 ArrayType m5 = ArrayType::Random(131, 131);
Rasmus Munk Larsenb47c7772020-04-27 18:55:15 -0700739 ArrayType m6 = m5;
740 m6.transposeInPlace();
741 VERIFY_IS_APPROX(m6, m5.transpose());
Jitse Niesen837db082011-05-09 10:17:41 +0100742}
743
Abraham Bachrach039408c2012-01-11 11:00:30 -0500744template<typename ArrayType> void min_max(const ArrayType& m)
745{
Abraham Bachrach039408c2012-01-11 11:00:30 -0500746 typedef typename ArrayType::Scalar Scalar;
747
748 Index rows = m.rows();
749 Index cols = m.cols();
750
751 ArrayType m1 = ArrayType::Random(rows, cols);
752
753 // min/max with array
754 Scalar maxM1 = m1.maxCoeff();
755 Scalar minM1 = m1.minCoeff();
756
757 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1)));
758 VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1)));
759
760 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1)));
761 VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1)));
762
763 // min/max with scalar input
764 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1));
765 VERIFY_IS_APPROX(m1, (m1.min)( maxM1));
766
767 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1));
768 VERIFY_IS_APPROX(m1, (m1.max)( minM1));
769
Rasmus Munk Larsen841c8982021-02-24 17:49:20 -0800770
771 // min/max with various NaN propagation options.
772 if (m1.size() > 1 && !NumTraits<Scalar>::IsInteger) {
Antonio Sanchez8dfe1022021-03-16 20:12:46 -0700773 m1(0,0) = NumTraits<Scalar>::quiet_NaN();
Rasmus Munk Larsen841c8982021-02-24 17:49:20 -0800774 maxM1 = m1.template maxCoeff<PropagateNaN>();
775 minM1 = m1.template minCoeff<PropagateNaN>();
776 VERIFY((numext::isnan)(maxM1));
777 VERIFY((numext::isnan)(minM1));
778
779 maxM1 = m1.template maxCoeff<PropagateNumbers>();
780 minM1 = m1.template minCoeff<PropagateNumbers>();
781 VERIFY(!(numext::isnan)(maxM1));
782 VERIFY(!(numext::isnan)(minM1));
783 }
Abraham Bachrach039408c2012-01-11 11:00:30 -0500784}
785
Antonio Sanchezfc9d3522021-08-17 20:04:48 -0700786template<int N>
787struct shift_left {
788 template<typename Scalar>
789 Scalar operator()(const Scalar& v) const {
790 return v << N;
791 }
792};
793
794template<int N>
795struct arithmetic_shift_right {
796 template<typename Scalar>
797 Scalar operator()(const Scalar& v) const {
798 return v >> N;
799 }
800};
801
802template<typename ArrayType> void array_integer(const ArrayType& m)
803{
804 Index rows = m.rows();
805 Index cols = m.cols();
806
807 ArrayType m1 = ArrayType::Random(rows, cols),
808 m2(rows, cols);
809
810 m2 = m1.template shiftLeft<2>();
811 VERIFY( (m2 == m1.unaryExpr(shift_left<2>())).all() );
812 m2 = m1.template shiftLeft<9>();
813 VERIFY( (m2 == m1.unaryExpr(shift_left<9>())).all() );
814
815 m2 = m1.template shiftRight<2>();
816 VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<2>())).all() );
817 m2 = m1.template shiftRight<9>();
818 VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<9>())).all() );
819}
820
Gael Guennebaude0f6d352018-09-20 18:07:32 +0200821EIGEN_DECLARE_TEST(array_cwise)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000822{
823 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100824 CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100825 CALL_SUBTEST_2( array(Array22f()) );
826 CALL_SUBTEST_3( array(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200827 CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
828 CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
829 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 +0100830 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 -0700831 CALL_SUBTEST_6( array_integer(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
832 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 +0000833 }
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100834 for(int i = 0; i < g_repeat; i++) {
835 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100836 CALL_SUBTEST_2( comparisons(Array22f()) );
837 CALL_SUBTEST_3( comparisons(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200838 CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
839 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 +0100840 }
841 for(int i = 0; i < g_repeat; i++) {
Abraham Bachrach039408c2012-01-11 11:00:30 -0500842 CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) );
843 CALL_SUBTEST_2( min_max(Array22f()) );
844 CALL_SUBTEST_3( min_max(Array44d()) );
845 CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
846 CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
847 }
848 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100849 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100850 CALL_SUBTEST_2( array_real(Array22f()) );
851 CALL_SUBTEST_3( array_real(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200852 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 -0800853 CALL_SUBTEST_7( array_real(Array<Eigen::half, 32, 32>()) );
854 CALL_SUBTEST_8( array_real(Array<Eigen::bfloat16, 32, 32>()) );
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100855 }
Jitse Niesen837db082011-05-09 10:17:41 +0100856 for(int i = 0; i < g_repeat; i++) {
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200857 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 +0100858 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400859
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000860 for(int i = 0; i < g_repeat; i++) {
861 CALL_SUBTEST_6( int_pow_test() );
862 CALL_SUBTEST_7( mixed_pow_test() );
863 }
864
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100865 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
866 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
867 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
Gael Guennebaud64fcfd32016-06-14 11:26:57 +0200868 typedef CwiseUnaryOp<internal::scalar_abs_op<double>, ArrayXd > Xpr;
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100869 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
870 ArrayBase<Xpr>
871 >::value));
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000872}