blob: b44837fd2f15b5aa54c17983ff658b9d3e46cc6e [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
Rasmus Munk Larsendceb7792022-09-12 15:51:27 -0700147template <typename Base, typename Exponent, bool ExpIsInteger = NumTraits<Exponent>::IsInteger>
148struct ref_pow {
149 static Base run(Base base, Exponent exponent) {
150 EIGEN_USING_STD(pow);
151 return pow(base, static_cast<Base>(exponent));
152 }
153};
154
155template <typename Base, typename Exponent>
156struct ref_pow<Base, Exponent, true> {
157 static Base run(Base base, Exponent exponent) {
158 EIGEN_USING_STD(pow);
159 return pow(base, exponent);
160 }
161};
162
Charles Schlossere5af9f82022-08-29 19:23:54 +0000163template <typename Base, typename Exponent>
164void test_exponent(Exponent exponent) {
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000165 const Base max_abs_bases = static_cast<Base>(10000);
166 // avoid integer overflow in Base type
167 Base threshold = calc_overflow_threshold<Base, Exponent>(numext::abs(exponent));
168 // avoid numbers that can't be verified with std::pow
169 double double_threshold = calc_overflow_threshold<double, Exponent>(numext::abs(exponent));
170 // use the lesser of these two thresholds
171 Base testing_threshold =
172 static_cast<double>(threshold) < double_threshold ? threshold : static_cast<Base>(double_threshold);
173 // test both vectorized and non-vectorized code paths
174 const Index array_size = 2 * internal::packet_traits<Base>::size + 1;
175
176 Base max_base = numext::mini(testing_threshold, max_abs_bases);
177 Base min_base = NumTraits<Base>::IsSigned ? -max_base : Base(0);
178
179 ArrayX<Base> x(array_size), y(array_size);
180 bool all_pass = true;
181 for (Base base = min_base; base <= max_base; base++) {
182 if (exponent < 0 && base == 0) continue;
183 x.setConstant(base);
184 y = x.pow(exponent);
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000185 for (Base a : y) {
Rasmus Munk Larsendceb7792022-09-12 15:51:27 -0700186 Base e = ref_pow<Base, Exponent>::run(base, exponent);
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000187 bool pass = (a == e);
188 if (!NumTraits<Base>::IsInteger) {
189 pass = pass || (((numext::isfinite)(e) && internal::isApprox(a, e)) ||
190 ((numext::isnan)(a) && (numext::isnan)(e)));
191 }
192 all_pass &= pass;
193 if (!pass) {
194 std::cout << "pow(" << base << "," << exponent << ") = " << a << " != " << e << std::endl;
195 }
Charles Schlossere5af9f82022-08-29 19:23:54 +0000196 }
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000197 }
198 VERIFY(all_pass);
Charles Schlossere5af9f82022-08-29 19:23:54 +0000199}
Charles Schlossere5af9f82022-08-29 19:23:54 +0000200
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000201template <typename Base, typename Exponent>
202void unary_pow_test() {
203 Exponent max_exponent = static_cast<Exponent>(NumTraits<Base>::digits());
204 Exponent min_exponent = static_cast<Exponent>(NumTraits<Exponent>::IsSigned ? -max_exponent : 0);
205
206 for (Exponent exponent = min_exponent; exponent < max_exponent; ++exponent) {
207 test_exponent<Base, Exponent>(exponent);
208 }
209};
210
211void mixed_pow_test() {
212 // The following cases will test promoting a smaller exponent type
213 // to a wider base type.
214 unary_pow_test<double, int>();
215 unary_pow_test<double, float>();
216 unary_pow_test<float, half>();
217 unary_pow_test<double, half>();
218 unary_pow_test<float, bfloat16>();
219 unary_pow_test<double, bfloat16>();
220
221 // Although in the following cases the exponent cannot be represented exactly
222 // in the base type, we do not perform a conversion, but implement
223 // the operation using repeated squaring.
224 unary_pow_test<float, int>();
225 unary_pow_test<double, long long>();
226
227 // The following cases will test promoting a wider exponent type
228 // to a narrower base type. This should compile but generate a
229 // deprecation warning:
230 unary_pow_test<float, double>();
231}
232
233void int_pow_test() {
234 unary_pow_test<int, int>();
235 unary_pow_test<unsigned int, unsigned int>();
236 unary_pow_test<long long, long long>();
237 unary_pow_test<unsigned long long, unsigned long long>();
238
239 // Although in the following cases the exponent cannot be represented exactly
240 // in the base type, we do not perform a conversion, but implement the
241 // operation using repeated squaring.
242 unary_pow_test<long long, int>();
243 unary_pow_test<int, unsigned int>();
244 unary_pow_test<unsigned int, int>();
245 unary_pow_test<long long, unsigned long long>();
246 unary_pow_test<unsigned long long, long long>();
247 unary_pow_test<long long, int>();
Charles Schlossere5af9f82022-08-29 19:23:54 +0000248}
249
Benoit Jacobe2775862010-04-28 18:51:38 -0400250template<typename ArrayType> void array(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000251{
Benoit Jacobe2775862010-04-28 18:51:38 -0400252 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200253 typedef typename ArrayType::RealScalar RealScalar;
Benoit Jacobe2775862010-04-28 18:51:38 -0400254 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
255 typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000256
Hauke Heibelf1679c72010-06-20 17:37:56 +0200257 Index rows = m.rows();
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800258 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000259
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000260 ArrayType m1 = ArrayType::Random(rows, cols);
261 if (NumTraits<RealScalar>::IsInteger && NumTraits<RealScalar>::IsSigned
262 && !NumTraits<Scalar>::IsComplex) {
263 // Here we cap the size of the values in m1 such that pow(3)/cube()
264 // doesn't overflow and result in undefined behavior. Notice that because
265 // pow(int, int) promotes its inputs and output to double (according to
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000266 // the C++ standard), we have to make sure that the result fits in 53 bits
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000267 // for int64,
268 RealScalar max_val =
269 numext::mini(RealScalar(std::cbrt(NumTraits<RealScalar>::highest())),
270 RealScalar(std::cbrt(1LL << 53)))/2;
271 m1.array() = (m1.abs().array() <= max_val).select(m1, Scalar(max_val));
272 }
273 ArrayType m2 = ArrayType::Random(rows, cols),
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000274 m3(rows, cols);
Christoph Hertzberg3be9f5c2015-04-16 13:25:20 +0200275 ArrayType m4 = m1; // copy constructor
276 VERIFY_IS_APPROX(m1, m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000277
Gael Guennebaud627595a2009-06-10 11:20:30 +0200278 ColVectorType cv1 = ColVectorType::Random(rows);
279 RowVectorType rv1 = RowVectorType::Random(cols);
280
Benoit Jacob47160402010-10-25 10:15:22 -0400281 Scalar s1 = internal::random<Scalar>(),
Hauke Heibel8cb3e362012-03-02 16:27:27 +0100282 s2 = internal::random<Scalar>();
283
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000284 // scalar addition
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100285 VERIFY_IS_APPROX(m1 + s1, s1 + m1);
Benoit Jacobe2775862010-04-28 18:51:38 -0400286 VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100287 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
Benoit Jacobe2775862010-04-28 18:51:38 -0400288 VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
289 VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
290 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000291 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100292 m3 += s2;
293 VERIFY_IS_APPROX(m3, m1 + s2);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000294 m3 = m1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100295 m3 -= s1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800296 VERIFY_IS_APPROX(m3, m1 - s1);
297
Hauke Heibelf578dc72010-12-16 17:34:13 +0100298 // scalar operators via Maps
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000299 m3 = m1; m4 = m1;
300 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
301 VERIFY_IS_APPROX(m4, m3 - m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800302
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000303 m3 = m1; m4 = m1;
304 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
305 VERIFY_IS_APPROX(m4, m3 + m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800306
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000307 m3 = m1; m4 = m1;
308 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
309 VERIFY_IS_APPROX(m4, m3 * m2);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800310
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000311 m3 = m1; m4 = m1;
Hauke Heibelf578dc72010-12-16 17:34:13 +0100312 m2 = ArrayType::Random(rows,cols);
313 m2 = (m2==0).select(1,m2);
Rasmus Munk Larsen98e51c92022-08-26 16:19:03 +0000314 ArrayType::Map(m4.data(), m4.rows(), m4.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
315 VERIFY_IS_APPROX(m4, m3 / m2);
Gael Guennebaud05ad0832008-07-19 13:03:23 +0000316
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000317 // reductions
Gael Guennebaud42af5872013-02-23 22:58:14 +0100318 VERIFY_IS_APPROX(m1.abs().colwise().sum().sum(), m1.abs().sum());
319 VERIFY_IS_APPROX(m1.abs().rowwise().sum().sum(), m1.abs().sum());
320 using std::abs;
321 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.colwise().sum().sum() - m1.sum()), m1.abs().sum());
322 VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.rowwise().sum().sum() - m1.sum()), m1.abs().sum());
Gael Guennebaud28e139a2013-02-23 23:06:45 +0100323 if (!internal::isMuchSmallerThan(abs(m1.sum() - (m1+m2).sum()), m1.abs().sum(), test_precision<Scalar>()))
Hauke Heibel83e3c452010-12-16 18:53:02 +0100324 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
Gael Guennebaud66e99ab2016-06-06 15:11:41 +0200325 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
Gael Guennebaud627595a2009-06-10 11:20:30 +0200326
327 // vector-wise ops
328 m3 = m1;
329 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
330 m3 = m1;
331 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
332 m3 = m1;
333 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
334 m3 = m1;
335 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800336
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200337 // Conversion from scalar
338 VERIFY_IS_APPROX((m3 = s1), ArrayType::Constant(rows,cols,s1));
339 VERIFY_IS_APPROX((m3 = 1), ArrayType::Constant(rows,cols,1));
340 VERIFY_IS_APPROX((m3.topLeftCorner(rows,cols) = 1), ArrayType::Constant(rows,cols,1));
341 typedef Array<Scalar,
342 ArrayType::RowsAtCompileTime==Dynamic?2:ArrayType::RowsAtCompileTime,
343 ArrayType::ColsAtCompileTime==Dynamic?2:ArrayType::ColsAtCompileTime,
344 ArrayType::Options> FixedArrayType;
Gael Guennebaud543529d2019-01-22 15:30:50 +0100345 {
346 FixedArrayType f1(s1);
347 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
348 FixedArrayType f2(numext::real(s1));
349 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
350 FixedArrayType f3((int)100*numext::real(s1));
351 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
352 f1.setRandom();
353 FixedArrayType f4(f1.data());
354 VERIFY_IS_APPROX(f4, f1);
355 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100356 {
357 FixedArrayType f1{s1};
358 VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
359 FixedArrayType f2{numext::real(s1)};
360 VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
361 FixedArrayType f3{(int)100*numext::real(s1)};
362 VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
363 f1.setRandom();
364 FixedArrayType f4{f1.data()};
365 VERIFY_IS_APPROX(f4, f1);
366 }
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800367
Gael Guennebaudd476cad2016-06-25 10:12:06 +0200368 // pow
369 VERIFY_IS_APPROX(m1.pow(2), m1.square());
370 VERIFY_IS_APPROX(pow(m1,2), m1.square());
371 VERIFY_IS_APPROX(m1.pow(3), m1.cube());
372 VERIFY_IS_APPROX(pow(m1,3), m1.cube());
373 VERIFY_IS_APPROX((-m1).pow(3), -m1.cube());
374 VERIFY_IS_APPROX(pow(2*m1,3), 8*m1.cube());
375 ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2));
376 VERIFY_IS_APPROX(Eigen::pow(m1,exponents), m1.square());
377 VERIFY_IS_APPROX(m1.pow(exponents), m1.square());
378 VERIFY_IS_APPROX(Eigen::pow(2*m1,exponents), 4*m1.square());
379 VERIFY_IS_APPROX((2*m1).pow(exponents), 4*m1.square());
380 VERIFY_IS_APPROX(Eigen::pow(m1,2*exponents), m1.square().square());
381 VERIFY_IS_APPROX(m1.pow(2*exponents), m1.square().square());
Gael Guennebaude61cee72016-07-04 11:49:03 +0200382 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 +0200383
Gael Guennebaud0a18eec2014-09-19 13:25:28 +0200384 // Check possible conflicts with 1D ctor
385 typedef Array<Scalar, Dynamic, 1> OneDArrayType;
Gael Guennebaud543529d2019-01-22 15:30:50 +0100386 {
387 OneDArrayType o1(rows);
388 VERIFY(o1.size()==rows);
389 OneDArrayType o2(static_cast<int>(rows));
390 VERIFY(o2.size()==rows);
391 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100392 {
393 OneDArrayType o1{rows};
394 VERIFY(o1.size()==rows);
395 OneDArrayType o4{int(rows)};
396 VERIFY(o4.size()==rows);
397 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100398 // Check possible conflicts with 2D ctor
399 typedef Array<Scalar, Dynamic, Dynamic> TwoDArrayType;
400 typedef Array<Scalar, 2, 1> ArrayType2;
401 {
402 TwoDArrayType o1(rows,cols);
403 VERIFY(o1.rows()==rows);
404 VERIFY(o1.cols()==cols);
405 TwoDArrayType o2(static_cast<int>(rows),static_cast<int>(cols));
406 VERIFY(o2.rows()==rows);
407 VERIFY(o2.cols()==cols);
408
409 ArrayType2 o3(rows,cols);
410 VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
411 ArrayType2 o4(static_cast<int>(rows),static_cast<int>(cols));
412 VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
413 }
Gael Guennebaud543529d2019-01-22 15:30:50 +0100414 {
415 TwoDArrayType o1{rows,cols};
416 VERIFY(o1.rows()==rows);
417 VERIFY(o1.cols()==cols);
418 TwoDArrayType o2{int(rows),int(cols)};
419 VERIFY(o2.rows()==rows);
420 VERIFY(o2.cols()==cols);
421
422 ArrayType2 o3{rows,cols};
423 VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols));
424 ArrayType2 o4{int(rows),int(cols)};
425 VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols));
426 }
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000427}
428
Benoit Jacobe2775862010-04-28 18:51:38 -0400429template<typename ArrayType> void comparisons(const ArrayType& m)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000430{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100431 using std::abs;
Benoit Jacobe2775862010-04-28 18:51:38 -0400432 typedef typename ArrayType::Scalar Scalar;
Benoit Jacob874ff5a2009-01-26 17:56:04 +0000433 typedef typename NumTraits<Scalar>::Real RealScalar;
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000434
Hauke Heibelf1679c72010-06-20 17:37:56 +0200435 Index rows = m.rows();
436 Index cols = m.cols();
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000437
Benoit Jacob47160402010-10-25 10:15:22 -0400438 Index r = internal::random<Index>(0, rows-1),
439 c = internal::random<Index>(0, cols-1);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000440
Benoit Jacobe2775862010-04-28 18:51:38 -0400441 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200442 m2 = ArrayType::Random(rows, cols),
443 m3(rows, cols),
444 m4 = m1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800445
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200446 m4 = (m4.abs()==Scalar(0)).select(1,m4);
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000447
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100448 VERIFY(((m1 + Scalar(1)) > m1).all());
449 VERIFY(((m1 - Scalar(1)) < m1).all());
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000450 if (rows*cols>1)
451 {
452 m3 = m1;
453 m3(r,c) += 1;
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100454 VERIFY(! (m1 < m3).all() );
455 VERIFY(! (m1 > m3).all() );
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000456 }
Christoph Hertzberg494fa992015-05-07 17:28:40 +0200457 VERIFY(!(m1 > m2 && m1 < m2).any());
458 VERIFY((m1 <= m2 || m1 >= m2).all());
Gael Guennebaud627595a2009-06-10 11:20:30 +0200459
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200460 // comparisons array to scalar
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100461 VERIFY( (m1 != (m1(r,c)+1) ).any() );
Christoph Hertzberg36052c42013-10-17 15:37:29 +0200462 VERIFY( (m1 > (m1(r,c)-1) ).any() );
463 VERIFY( (m1 < (m1(r,c)+1) ).any() );
464 VERIFY( (m1 == m1(r,c) ).any() );
465
466 // comparisons scalar to array
467 VERIFY( ( (m1(r,c)+1) != m1).any() );
468 VERIFY( ( (m1(r,c)-1) < m1).any() );
469 VERIFY( ( (m1(r,c)+1) > m1).any() );
470 VERIFY( ( m1(r,c) == m1).any() );
Gael Guennebaud627595a2009-06-10 11:20:30 +0200471
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000472 // test Select
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100473 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
474 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
Gael Guennebaud9f795582009-12-16 19:18:40 +0100475 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000476 for (int j=0; j<cols; ++j)
477 for (int i=0; i<rows; ++i)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100478 m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
Benoit Jacobe2775862010-04-28 18:51:38 -0400479 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
480 .select(ArrayType::Zero(rows,cols),m1), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000481 // shorter versions:
Benoit Jacobe2775862010-04-28 18:51:38 -0400482 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000483 .select(0,m1), m3);
Benoit Jacobe2775862010-04-28 18:51:38 -0400484 VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
Gael Guennebaud59dc1da2008-09-03 17:16:28 +0000485 .select(m1,0), m3);
Gael Guennebaude14aa8c2008-09-03 17:56:06 +0000486 // even shorter version:
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100487 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
Gael Guennebaud627595a2009-06-10 11:20:30 +0200488
Gael Guennebaud56c7e162009-01-24 15:22:44 +0000489 // count
Gael Guennebaudc70d5422010-01-18 22:54:20 +0100490 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
Benoit Jacobaaaade42010-05-30 16:00:58 -0400491
Gael Guennebaud35c11582011-05-31 22:17:34 +0200492 // and/or
493 VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0);
494 VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols);
495 RealScalar a = m1.abs().mean();
496 VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count());
497
Gael Guennebaud12e1ebb2018-07-12 17:16:40 +0200498 typedef Array<Index, Dynamic, 1> ArrayOfIndices;
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200499
Gael Guennebaud9f795582009-12-16 19:18:40 +0100500 // TODO allows colwise/rowwise for array
Benoit Jacobaaaade42010-05-30 16:00:58 -0400501 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
502 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
Benoit Jacobe8009992008-11-03 22:47:00 +0000503}
504
Benoit Jacobe2775862010-04-28 18:51:38 -0400505template<typename ArrayType> void array_real(const ArrayType& m)
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100506{
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100507 using std::abs;
Gael Guennebaud9b9177f2013-07-05 13:35:34 +0200508 using std::sqrt;
Benoit Jacobe2775862010-04-28 18:51:38 -0400509 typedef typename ArrayType::Scalar Scalar;
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100510 typedef typename NumTraits<Scalar>::Real RealScalar;
511
Hauke Heibelf1679c72010-06-20 17:37:56 +0200512 Index rows = m.rows();
513 Index cols = m.cols();
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100514
Benoit Jacobe2775862010-04-28 18:51:38 -0400515 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200516 m2 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200517 m3(rows, cols),
518 m4 = m1;
Benoit Steiner83149622015-12-10 13:13:45 -0800519
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800520 m4 = (m4.abs()==Scalar(0)).select(Scalar(1),m4);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100521
Gael Guennebaud9b1ad5e2012-03-09 12:08:06 +0100522 Scalar s1 = internal::random<Scalar>();
523
Deanna Hood41b717d2015-03-18 03:11:03 +1000524 // these tests are mostly to check possible compilation issues with free-functions.
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100525 VERIFY_IS_APPROX(m1.sin(), sin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100526 VERIFY_IS_APPROX(m1.cos(), cos(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000527 VERIFY_IS_APPROX(m1.tan(), tan(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100528 VERIFY_IS_APPROX(m1.asin(), asin(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100529 VERIFY_IS_APPROX(m1.acos(), acos(m1));
Jitse Niesende150b12014-06-19 15:12:33 +0100530 VERIFY_IS_APPROX(m1.atan(), atan(m1));
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000531 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
532 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
533 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Rasmus Munk Larsen1e1848f2022-09-28 20:46:49 +0000534 VERIFY_IS_APPROX(m1.atan2(m2), atan2(m1,m2));
535
Rasmus Munk Larsen28ba1b22019-01-11 17:45:37 -0800536#if EIGEN_HAS_CXX11_MATH
537 VERIFY_IS_APPROX(m1.tanh().atanh(), atanh(tanh(m1)));
538 VERIFY_IS_APPROX(m1.sinh().asinh(), asinh(sinh(m1)));
539 VERIFY_IS_APPROX(m1.cosh().acosh(), acosh(cosh(m1)));
540#endif
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700541 VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
Gael Guennebaud2f7e2612016-07-08 11:13:55 +0200542
Deanna Hooda5e49972015-03-11 08:56:42 +1000543 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000544 VERIFY_IS_APPROX(m1.round(), round(m1));
Ilya Tokar19876ce2019-12-16 16:00:35 -0500545 VERIFY_IS_APPROX(m1.rint(), rint(m1));
Deanna Hood31fdd672015-03-11 06:39:23 +1000546 VERIFY_IS_APPROX(m1.floor(), floor(m1));
547 VERIFY_IS_APPROX(m1.ceil(), ceil(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200548 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
549 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
550 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800551 VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
Deanna Hood41b717d2015-03-18 03:11:03 +1000552 VERIFY_IS_APPROX(m1.abs(), abs(m1));
553 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
554 VERIFY_IS_APPROX(m1.square(), square(m1));
555 VERIFY_IS_APPROX(m1.cube(), cube(m1));
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100556 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100557 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Rasmus Munk Larsen6964ae82020-07-07 01:54:04 +0000558 VERIFY((m1.sqrt().sign().isNaN() == (Eigen::isnan)(sign(sqrt(m1)))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000559
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800560 // avoid inf and NaNs so verification doesn't fail
561 m3 = m4.abs();
562 VERIFY_IS_APPROX(m3.sqrt(), sqrt(abs(m3)));
563 VERIFY_IS_APPROX(m3.rsqrt(), Scalar(1)/sqrt(abs(m3)));
564 VERIFY_IS_APPROX(rsqrt(m3), Scalar(1)/sqrt(abs(m3)));
Deanna Hood41b717d2015-03-18 03:11:03 +1000565 VERIFY_IS_APPROX(m3.log(), log(m3));
Gael Guennebaud89099b02016-06-01 17:00:08 +0200566 VERIFY_IS_APPROX(m3.log1p(), log1p(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000567 VERIFY_IS_APPROX(m3.log10(), log10(m3));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000568 VERIFY_IS_APPROX(m3.log2(), log2(m3));
Deanna Hood41b717d2015-03-18 03:11:03 +1000569
570
571 VERIFY((!(m1>m2) == (m1<=m2)).all());
572
573 VERIFY_IS_APPROX(sin(m1.asin()), m1);
574 VERIFY_IS_APPROX(cos(m1.acos()), m1);
575 VERIFY_IS_APPROX(tan(m1.atan()), m1);
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800576 VERIFY_IS_APPROX(sinh(m1), Scalar(0.5)*(exp(m1)-exp(-m1)));
577 VERIFY_IS_APPROX(cosh(m1), Scalar(0.5)*(exp(m1)+exp(-m1)));
578 VERIFY_IS_APPROX(tanh(m1), (Scalar(0.5)*(exp(m1)-exp(-m1)))/(Scalar(0.5)*(exp(m1)+exp(-m1))));
579 VERIFY_IS_APPROX(logistic(m1), (Scalar(1)/(Scalar(1)+exp(-m1))));
580 VERIFY_IS_APPROX(arg(m1), ((m1<Scalar(0)).template cast<Scalar>())*Scalar(std::acos(Scalar(-1))));
Deanna Hood41b717d2015-03-18 03:11:03 +1000581 VERIFY((round(m1) <= ceil(m1) && round(m1) >= floor(m1)).all());
Ilya Tokar19876ce2019-12-16 16:00:35 -0500582 VERIFY((rint(m1) <= ceil(m1) && rint(m1) >= floor(m1)).all());
583 VERIFY(((ceil(m1) - round(m1)) <= Scalar(0.5) || (round(m1) - floor(m1)) <= Scalar(0.5)).all());
584 VERIFY(((ceil(m1) - round(m1)) <= Scalar(1.0) && (round(m1) - floor(m1)) <= Scalar(1.0)).all());
585 VERIFY(((ceil(m1) - rint(m1)) <= Scalar(0.5) || (rint(m1) - floor(m1)) <= Scalar(0.5)).all());
586 VERIFY(((ceil(m1) - rint(m1)) <= Scalar(1.0) && (rint(m1) - floor(m1)) <= Scalar(1.0)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800587 VERIFY((Eigen::isnan)((m1*Scalar(0))/Scalar(0)).all());
588 VERIFY((Eigen::isinf)(m4/Scalar(0)).all());
589 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*Scalar(0)/Scalar(0))) && (!(Eigen::isfinite)(m4/Scalar(0)))).all());
590 VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
Deanna Hood41b717d2015-03-18 03:11:03 +1000591 VERIFY((abs(m1) == m1 || abs(m1) == -m1).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800592 VERIFY_IS_APPROX(m3, sqrt(abs2(m3)));
Joel Holdsworthd5c66572020-03-19 17:45:20 +0000593 VERIFY_IS_APPROX(m1.absolute_difference(m2), (m1 > m2).select(m1 - m2, m2 - m1));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100594 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
595 VERIFY_IS_APPROX( m1*m1.sign(),m1.abs());
596 VERIFY_IS_APPROX(m1.sign() * m1.abs(), m1);
Rasmus Munk Larsen1e1848f2022-09-28 20:46:49 +0000597
598 ArrayType tmp = m1.atan2(m2);
599 for (Index i = 0; i < tmp.size(); ++i) {
600 Scalar actual = tmp.array()(i);
601 Scalar expected = atan2(m1.array()(i), m2.array()(i));
602 VERIFY_IS_APPROX(actual, expected);
603 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400604
Gael Guennebaud62670c82013-06-10 23:40:56 +0200605 VERIFY_IS_APPROX(numext::abs2(numext::real(m1)) + numext::abs2(numext::imag(m1)), numext::abs2(m1));
Gael Guennebaud87427d22019-10-08 09:15:17 +0200606 VERIFY_IS_APPROX(numext::abs2(Eigen::real(m1)) + numext::abs2(Eigen::imag(m1)), numext::abs2(m1));
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400607 if(!NumTraits<Scalar>::IsComplex)
Gael Guennebaud62670c82013-06-10 23:40:56 +0200608 VERIFY_IS_APPROX(numext::real(m1), m1);
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200609
Jitse Niesen6fecb6f2014-02-24 14:10:17 +0000610 // shift argument of logarithm so that it is not zero
611 Scalar smallNumber = NumTraits<Scalar>::dummy_precision();
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800612 VERIFY_IS_APPROX((m3 + smallNumber).log() , log(abs(m3) + smallNumber));
613 VERIFY_IS_APPROX((m3 + smallNumber + Scalar(1)).log() , log1p(abs(m3) + smallNumber));
Gael Guennebaud4ebb8042010-06-02 09:45:57 +0200614
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100615 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
616 VERIFY_IS_APPROX(m1.exp(), exp(m1));
617 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
Gael Guennebaud575ac542010-06-19 23:17:07 +0200618
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800619 VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800620 VERIFY_IS_APPROX((m3 + smallNumber).exp() - Scalar(1), expm1(abs(m3) + smallNumber));
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800621
Gael Guennebaud575ac542010-06-19 23:17:07 +0200622 VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100623 VERIFY_IS_APPROX(pow(m3,RealScalar(0.5)), m3.sqrt());
Benoit Steinere25e3a02015-12-03 18:16:35 -0800624
625 VERIFY_IS_APPROX(m3.pow(RealScalar(-0.5)), m3.rsqrt());
626 VERIFY_IS_APPROX(pow(m3,RealScalar(-0.5)), m3.rsqrt());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800627
628 // Avoid inf and NaN.
629 m3 = (m1.square()<NumTraits<Scalar>::epsilon()).select(Scalar(1),m3);
630 VERIFY_IS_APPROX(m3.pow(RealScalar(-2)), m3.square().inverse());
Rasmus Munk Larsencdd8fdc2021-01-18 13:25:16 +0000631 pow_test<Scalar>();
Benoit Steinere25e3a02015-12-03 18:16:35 -0800632
Antonio Sanchez4c42d5e2021-01-23 11:54:00 -0800633 VERIFY_IS_APPROX(log10(m3), log(m3)/numext::log(Scalar(10)));
634 VERIFY_IS_APPROX(log2(m3), log(m3)/numext::log(Scalar(2)));
Hauke Heibel81c13362012-03-07 08:58:42 +0100635
636 // scalar by array division
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100637 const RealScalar tiny = sqrt(std::numeric_limits<RealScalar>::epsilon());
Hauke Heibeldd9365e2012-03-09 14:04:13 +0100638 s1 += Scalar(tiny);
639 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
Rasmus Munk Larsen96dc37a2022-01-07 01:10:17 +0000640 VERIFY_IS_CWISE_APPROX(s1/m1, s1 * m1.inverse());
Eugene Brevdof7362772015-12-24 21:15:38 -0800641
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200642 // check inplace transpose
643 m3 = m1;
644 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000645 VERIFY_IS_APPROX(m3, m1.transpose());
Gael Guennebaudc695cbf2013-06-24 13:33:44 +0200646 m3.transposeInPlace();
Deanna Hood41b717d2015-03-18 03:11:03 +1000647 VERIFY_IS_APPROX(m3, m1);
Gael Guennebaud0ce5bc02010-01-27 23:23:59 +0100648}
649
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000650
Jitse Niesen837db082011-05-09 10:17:41 +0100651template<typename ArrayType> void array_complex(const ArrayType& m)
652{
Deanna Hood41b717d2015-03-18 03:11:03 +1000653 typedef typename ArrayType::Scalar Scalar;
654 typedef typename NumTraits<Scalar>::Real RealScalar;
Jitse Niesen837db082011-05-09 10:17:41 +0100655
656 Index rows = m.rows();
657 Index cols = m.cols();
658
659 ArrayType m1 = ArrayType::Random(rows, cols),
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200660 m2(rows, cols),
661 m4 = m1;
Srinivas Vasudevan218764e2016-12-02 14:13:01 -0800662
Gael Guennebaud9fc1c922015-06-22 16:48:27 +0200663 m4.real() = (m4.real().abs()==RealScalar(0)).select(RealScalar(1),m4.real());
664 m4.imag() = (m4.imag().abs()==RealScalar(0)).select(RealScalar(1),m4.imag());
Jitse Niesen837db082011-05-09 10:17:41 +0100665
Deanna Hood41b717d2015-03-18 03:11:03 +1000666 Array<RealScalar, -1, -1> m3(rows, cols);
667
Jitse Niesen837db082011-05-09 10:17:41 +0100668 for (Index i = 0; i < m.rows(); ++i)
669 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebauda76fbbf2012-11-06 15:25:50 +0100670 m2(i,j) = sqrt(m1(i,j));
Jitse Niesen837db082011-05-09 10:17:41 +0100671
Deanna Hood41b717d2015-03-18 03:11:03 +1000672 // these tests are mostly to check possible compilation issues with free-functions.
Deanna Hoodf89fcef2015-03-11 13:13:30 +1000673 VERIFY_IS_APPROX(m1.sin(), sin(m1));
674 VERIFY_IS_APPROX(m1.cos(), cos(m1));
675 VERIFY_IS_APPROX(m1.tan(), tan(m1));
676 VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
677 VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
678 VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
Rasmus Munk Larsend6e283b2018-08-13 11:14:50 -0700679 VERIFY_IS_APPROX(m1.logistic(), logistic(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000680 VERIFY_IS_APPROX(m1.arg(), arg(m1));
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200681 VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
682 VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
683 VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800684 VERIFY_IS_APPROX(m4.inverse(), inverse(m4));
Deanna Hood41b717d2015-03-18 03:11:03 +1000685 VERIFY_IS_APPROX(m1.log(), log(m1));
686 VERIFY_IS_APPROX(m1.log10(), log10(m1));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000687 VERIFY_IS_APPROX(m1.log2(), log2(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000688 VERIFY_IS_APPROX(m1.abs(), abs(m1));
689 VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
690 VERIFY_IS_APPROX(m1.sqrt(), sqrt(m1));
691 VERIFY_IS_APPROX(m1.square(), square(m1));
692 VERIFY_IS_APPROX(m1.cube(), cube(m1));
693 VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100694 VERIFY_IS_APPROX(m1.sign(), sign(m1));
Deanna Hood41b717d2015-03-18 03:11:03 +1000695
Deanna Hood41b717d2015-03-18 03:11:03 +1000696 VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
697 VERIFY_IS_APPROX(m1.exp(), exp(m1));
698 VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
699
Srinivas Vasudevan88062b72019-12-12 01:56:54 +0000700 VERIFY_IS_APPROX(m1.expm1(), expm1(m1));
701 VERIFY_IS_APPROX(expm1(m1), exp(m1) - 1.);
702 // Check for larger magnitude complex numbers that expm1 matches exp - 1.
703 VERIFY_IS_APPROX(expm1(10. * m1), exp(10. * m1) - 1.);
704
Deanna Hood41b717d2015-03-18 03:11:03 +1000705 VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
706 VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
707 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 -0700708 VERIFY_IS_APPROX(logistic(m1), (1.0/(1.0 + exp(-m1))));
Deanna Hood41b717d2015-03-18 03:11:03 +1000709
710 for (Index i = 0; i < m.rows(); ++i)
711 for (Index j = 0; j < m.cols(); ++j)
Gael Guennebaud87427d22019-10-08 09:15:17 +0200712 m3(i,j) = std::atan2(m1(i,j).imag(), m1(i,j).real());
Deanna Hood41b717d2015-03-18 03:11:03 +1000713 VERIFY_IS_APPROX(arg(m1), m3);
714
715 std::complex<RealScalar> zero(0.0,0.0);
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200716 VERIFY((Eigen::isnan)(m1*zero/zero).all());
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100717#if EIGEN_COMP_MSVC
718 // msvc complex division is not robust
719 VERIFY((Eigen::isinf)(m4/RealScalar(0)).all());
720#else
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200721#if EIGEN_COMP_CLANG
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100722 // clang's complex division is notoriously broken too
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200723 if((numext::isinf)(m4(0,0)/RealScalar(0))) {
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200724#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100725 VERIFY((Eigen::isinf)(m4/zero).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200726#if EIGEN_COMP_CLANG
727 }
728 else
729 {
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200730 VERIFY((Eigen::isinf)(m4.real()/zero.real()).all());
Gael Guennebaud40f326e2015-06-17 15:33:09 +0200731 }
732#endif
Gael Guennebaud4a985e72015-11-20 14:52:08 +0100733#endif // MSVC
734
Christoph Hertzberg61e09772015-08-14 17:32:34 +0200735 VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*zero/zero)) && (!(Eigen::isfinite)(m1/zero))).all());
Deanna Hood41b717d2015-03-18 03:11:03 +1000736
Antonio Sanchezf0e46ed2021-01-22 11:10:54 -0800737 VERIFY_IS_APPROX(inverse(inverse(m4)),m4);
Deanna Hood41b717d2015-03-18 03:11:03 +1000738 VERIFY_IS_APPROX(conj(m1.conjugate()), m1);
Gael Guennebaud87427d22019-10-08 09:15:17 +0200739 VERIFY_IS_APPROX(abs(m1), sqrt(square(m1.real())+square(m1.imag())));
Deanna Hood41b717d2015-03-18 03:11:03 +1000740 VERIFY_IS_APPROX(abs(m1), sqrt(abs2(m1)));
741 VERIFY_IS_APPROX(log10(m1), log(m1)/log(10));
Rasmus Munk Larsenf9fac1d2020-12-04 21:45:09 +0000742 VERIFY_IS_APPROX(log2(m1), log(m1)/log(2));
Deanna Hood41b717d2015-03-18 03:11:03 +1000743
Gael Guennebaud3f32f5e2015-11-27 16:27:53 +0100744 VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
745 VERIFY_IS_APPROX( m1.sign() * m1.abs(), m1);
746
Deanna Hood41b717d2015-03-18 03:11:03 +1000747 // scalar by array division
Eugene Brevdof7362772015-12-24 21:15:38 -0800748 Scalar s1 = internal::random<Scalar>();
Benoit Steiner2d74ef92016-05-17 07:20:11 -0700749 const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon());
Deanna Hood41b717d2015-03-18 03:11:03 +1000750 s1 += Scalar(tiny);
751 m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
752 VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
753
754 // check inplace transpose
755 m2 = m1;
756 m2.transposeInPlace();
757 VERIFY_IS_APPROX(m2, m1.transpose());
758 m2.transposeInPlace();
759 VERIFY_IS_APPROX(m2, m1);
Rasmus Munk Larsenb47c7772020-04-27 18:55:15 -0700760 // Check vectorized inplace transpose.
Rasmus Munk Larsen74ec8e62020-05-07 17:29:56 +0000761 ArrayType m5 = ArrayType::Random(131, 131);
Rasmus Munk Larsenb47c7772020-04-27 18:55:15 -0700762 ArrayType m6 = m5;
763 m6.transposeInPlace();
764 VERIFY_IS_APPROX(m6, m5.transpose());
Jitse Niesen837db082011-05-09 10:17:41 +0100765}
766
Abraham Bachrach039408c2012-01-11 11:00:30 -0500767template<typename ArrayType> void min_max(const ArrayType& m)
768{
Abraham Bachrach039408c2012-01-11 11:00:30 -0500769 typedef typename ArrayType::Scalar Scalar;
770
771 Index rows = m.rows();
772 Index cols = m.cols();
773
774 ArrayType m1 = ArrayType::Random(rows, cols);
775
776 // min/max with array
777 Scalar maxM1 = m1.maxCoeff();
778 Scalar minM1 = m1.minCoeff();
779
780 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1)));
781 VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1)));
782
783 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1)));
784 VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1)));
785
786 // min/max with scalar input
787 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1));
788 VERIFY_IS_APPROX(m1, (m1.min)( maxM1));
789
790 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1));
791 VERIFY_IS_APPROX(m1, (m1.max)( minM1));
792
Rasmus Munk Larsen841c8982021-02-24 17:49:20 -0800793
794 // min/max with various NaN propagation options.
795 if (m1.size() > 1 && !NumTraits<Scalar>::IsInteger) {
Antonio Sanchez8dfe1022021-03-16 20:12:46 -0700796 m1(0,0) = NumTraits<Scalar>::quiet_NaN();
Rasmus Munk Larsen841c8982021-02-24 17:49:20 -0800797 maxM1 = m1.template maxCoeff<PropagateNaN>();
798 minM1 = m1.template minCoeff<PropagateNaN>();
799 VERIFY((numext::isnan)(maxM1));
800 VERIFY((numext::isnan)(minM1));
801
802 maxM1 = m1.template maxCoeff<PropagateNumbers>();
803 minM1 = m1.template minCoeff<PropagateNumbers>();
804 VERIFY(!(numext::isnan)(maxM1));
805 VERIFY(!(numext::isnan)(minM1));
806 }
Abraham Bachrach039408c2012-01-11 11:00:30 -0500807}
808
Antonio Sanchezfc9d3522021-08-17 20:04:48 -0700809template<int N>
810struct shift_left {
811 template<typename Scalar>
812 Scalar operator()(const Scalar& v) const {
813 return v << N;
814 }
815};
816
817template<int N>
818struct arithmetic_shift_right {
819 template<typename Scalar>
820 Scalar operator()(const Scalar& v) const {
821 return v >> N;
822 }
823};
824
825template<typename ArrayType> void array_integer(const ArrayType& m)
826{
827 Index rows = m.rows();
828 Index cols = m.cols();
829
830 ArrayType m1 = ArrayType::Random(rows, cols),
831 m2(rows, cols);
832
833 m2 = m1.template shiftLeft<2>();
834 VERIFY( (m2 == m1.unaryExpr(shift_left<2>())).all() );
835 m2 = m1.template shiftLeft<9>();
836 VERIFY( (m2 == m1.unaryExpr(shift_left<9>())).all() );
837
838 m2 = m1.template shiftRight<2>();
839 VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<2>())).all() );
840 m2 = m1.template shiftRight<9>();
841 VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<9>())).all() );
842}
843
Gael Guennebaude0f6d352018-09-20 18:07:32 +0200844EIGEN_DECLARE_TEST(array_cwise)
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000845{
846 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100847 CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100848 CALL_SUBTEST_2( array(Array22f()) );
849 CALL_SUBTEST_3( array(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200850 CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
851 CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
852 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 +0100853 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 -0700854 CALL_SUBTEST_6( array_integer(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
855 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 +0000856 }
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100857 for(int i = 0; i < g_repeat; i++) {
858 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100859 CALL_SUBTEST_2( comparisons(Array22f()) );
860 CALL_SUBTEST_3( comparisons(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200861 CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
862 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 +0100863 }
864 for(int i = 0; i < g_repeat; i++) {
Abraham Bachrach039408c2012-01-11 11:00:30 -0500865 CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) );
866 CALL_SUBTEST_2( min_max(Array22f()) );
867 CALL_SUBTEST_3( min_max(Array44d()) );
868 CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
869 CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
870 }
871 for(int i = 0; i < g_repeat; i++) {
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100872 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
Hauke Heibelb31e1242010-12-16 19:07:23 +0100873 CALL_SUBTEST_2( array_real(Array22f()) );
874 CALL_SUBTEST_3( array_real(Array44d()) );
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200875 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 -0800876 CALL_SUBTEST_7( array_real(Array<Eigen::half, 32, 32>()) );
877 CALL_SUBTEST_8( array_real(Array<Eigen::bfloat16, 32, 32>()) );
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100878 }
Jitse Niesen837db082011-05-09 10:17:41 +0100879 for(int i = 0; i < g_repeat; i++) {
Gael Guennebauda8f66fe2011-07-12 14:41:00 +0200880 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 +0100881 }
Benoit Jacob5d63d2c2010-04-28 22:42:34 -0400882
Rasmus Munk Larsenafc014f2022-09-12 21:55:30 +0000883 for(int i = 0; i < g_repeat; i++) {
884 CALL_SUBTEST_6( int_pow_test() );
885 CALL_SUBTEST_7( mixed_pow_test() );
886 }
887
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100888 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
889 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
890 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
Gael Guennebaud64fcfd32016-06-14 11:26:57 +0200891 typedef CwiseUnaryOp<internal::scalar_abs_op<double>, ArrayXd > Xpr;
Hauke Heibel2d0dfe52010-12-16 17:36:10 +0100892 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
893 ArrayBase<Xpr>
894 >::value));
Gael Guennebaud8cef5412008-06-21 17:28:07 +0000895}