Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 1 | // This file is part of Eigen, a lightweight C++ template library |
Benoit Jacob | 6347b1d | 2009-05-22 20:25:33 +0200 | [diff] [blame] | 2 | // for linear algebra. |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 3 | // |
Gael Guennebaud | 28e64b0 | 2010-06-24 23:21:58 +0200 | [diff] [blame] | 4 | // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud@inria.fr> |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 5 | // |
Benoit Jacob | 69124cf | 2012-07-13 14:42:47 -0400 | [diff] [blame] | 6 | // This Source Code Form is subject to the terms of the Mozilla |
| 7 | // Public License v. 2.0. If a copy of the MPL was not distributed |
| 8 | // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. |
Gael Guennebaud | 19bc418 | 2013-02-25 01:17:44 +0100 | [diff] [blame] | 9 | |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 10 | #include "main.h" |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 11 | |
Rasmus Munk Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 12 | |
| 13 | // Test the corner cases of pow(x, y) for real types. |
| 14 | template<typename Scalar> |
| 15 | void pow_test() { |
| 16 | const Scalar zero = Scalar(0); |
Antonio Sanchez | 8dfe102 | 2021-03-16 20:12:46 -0700 | [diff] [blame] | 17 | const Scalar eps = Eigen::NumTraits<Scalar>::epsilon(); |
Rasmus Munk Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 18 | const Scalar one = Scalar(1); |
Rasmus Munk Larsen | 6e3b795 | 2021-02-05 16:58:49 -0800 | [diff] [blame] | 19 | const Scalar two = Scalar(2); |
| 20 | const Scalar three = Scalar(3); |
Rasmus Munk Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 21 | const Scalar sqrt_half = Scalar(std::sqrt(0.5)); |
| 22 | const Scalar sqrt2 = Scalar(std::sqrt(2)); |
Antonio Sanchez | 8dfe102 | 2021-03-16 20:12:46 -0700 | [diff] [blame] | 23 | const Scalar inf = Eigen::NumTraits<Scalar>::infinity(); |
| 24 | const Scalar nan = Eigen::NumTraits<Scalar>::quiet_NaN(); |
Rasmus Munk Larsen | be0574e | 2021-02-17 02:50:32 +0000 | [diff] [blame] | 25 | const Scalar denorm_min = std::numeric_limits<Scalar>::denorm_min(); |
Rasmus Munk Larsen | 6e3b795 | 2021-02-05 16:58:49 -0800 | [diff] [blame] | 26 | const Scalar min = (std::numeric_limits<Scalar>::min)(); |
| 27 | const Scalar max = (std::numeric_limits<Scalar>::max)(); |
Antonio Sanchez | 8dfe102 | 2021-03-16 20:12:46 -0700 | [diff] [blame] | 28 | const Scalar max_exp = (static_cast<Scalar>(int(Eigen::NumTraits<Scalar>::max_exponent())) * Scalar(EIGEN_LN2)) / eps; |
Rasmus Munk Larsen | be0574e | 2021-02-17 02:50:32 +0000 | [diff] [blame] | 29 | |
Rasmus Munk Larsen | 6e3b795 | 2021-02-05 16:58:49 -0800 | [diff] [blame] | 30 | const static Scalar abs_vals[] = {zero, |
Rasmus Munk Larsen | be0574e | 2021-02-17 02:50:32 +0000 | [diff] [blame] | 31 | denorm_min, |
| 32 | min, |
| 33 | eps, |
Rasmus Munk Larsen | 6e3b795 | 2021-02-05 16:58:49 -0800 | [diff] [blame] | 34 | sqrt_half, |
| 35 | one, |
| 36 | sqrt2, |
| 37 | two, |
| 38 | three, |
Rasmus Munk Larsen | be0574e | 2021-02-17 02:50:32 +0000 | [diff] [blame] | 39 | max_exp, |
Rasmus Munk Larsen | 6e3b795 | 2021-02-05 16:58:49 -0800 | [diff] [blame] | 40 | max, |
| 41 | inf, |
| 42 | nan}; |
Rasmus Munk Larsen | be0574e | 2021-02-17 02:50:32 +0000 | [diff] [blame] | 43 | const int abs_cases = 13; |
Rasmus Munk Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 44 | 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 Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 49 | 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 Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 61 | } |
| 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 Larsen | be0574e | 2021-02-17 02:50:32 +0000 | [diff] [blame] | 73 | Scalar e = static_cast<Scalar>(std::pow(x(i,j), y(i,j))); |
Rasmus Munk Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 74 | Scalar a = actual(i, j); |
Rasmus Munk Larsen | 7a87ed1 | 2022-08-08 18:48:36 +0000 | [diff] [blame] | 75 | 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 Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 78 | std::cout << "pow(" << x(i,j) << "," << y(i,j) << ") = " << a << " != " << e << std::endl; |
| 79 | } |
| 80 | } |
| 81 | } |
Charles Schlosser | 76a669f | 2022-08-16 21:32:36 +0000 | [diff] [blame] | 82 | |
| 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 Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 126 | VERIFY(all_pass); |
| 127 | } |
| 128 | |
Charles Schlosser | e5af9f8 | 2022-08-29 19:23:54 +0000 | [diff] [blame] | 129 | template <typename Scalar, typename ScalarExponent> |
| 130 | Scalar calc_overflow_threshold(const ScalarExponent exponent) { |
| 131 | EIGEN_USING_STD(exp2); |
| 132 | EIGEN_STATIC_ASSERT((NumTraits<Scalar>::digits() < 2 * NumTraits<double>::digits()), BASE_TYPE_IS_TOO_BIG); |
| 133 | |
| 134 | if (exponent < 2) |
| 135 | return NumTraits<Scalar>::highest(); |
| 136 | else { |
| 137 | const double max_exponent = static_cast<double>(NumTraits<Scalar>::digits()); |
| 138 | const double clamped_exponent = exponent < max_exponent ? static_cast<double>(exponent) : max_exponent; |
| 139 | const double threshold = exp2(max_exponent / clamped_exponent); |
| 140 | return static_cast<Scalar>(threshold); |
| 141 | } |
| 142 | } |
| 143 | |
| 144 | template <typename Base, typename Exponent> |
| 145 | void test_exponent(Exponent exponent) { |
| 146 | EIGEN_USING_STD(pow); |
| 147 | |
| 148 | const Base max_abs_bases = 10000; |
| 149 | // avoid integer overflow in Base type |
| 150 | Base threshold = calc_overflow_threshold<Base, Exponent>(numext::abs(exponent)); |
| 151 | // avoid numbers that can't be verified with std::pow |
| 152 | double double_threshold = calc_overflow_threshold<double, Exponent>(numext::abs(exponent)); |
| 153 | // use the lesser of these two thresholds |
| 154 | Base testing_threshold = threshold < double_threshold ? threshold : static_cast<Base>(double_threshold); |
Charles Schlosser | e5af9f8 | 2022-08-29 19:23:54 +0000 | [diff] [blame] | 155 | // test both vectorized and non-vectorized code paths |
| 156 | const Index array_size = 2 * internal::packet_traits<Base>::size + 1; |
| 157 | |
| 158 | Base max_base = numext::mini(testing_threshold, max_abs_bases); |
| 159 | Base min_base = NumTraits<Base>::IsSigned ? -max_base : 0; |
| 160 | |
| 161 | ArrayX<Base> x(array_size), y(array_size); |
| 162 | |
| 163 | bool all_pass = true; |
| 164 | |
| 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); |
| 169 | Base e = pow(base, exponent); |
| 170 | for (Base a : y) { |
| 171 | bool pass = a == e; |
| 172 | all_pass &= pass; |
| 173 | if (!pass) { |
| 174 | std::cout << "pow(" << base << "," << exponent << ") = " << a << " != " << e << std::endl; |
| 175 | } |
| 176 | } |
| 177 | } |
| 178 | |
| 179 | VERIFY(all_pass); |
| 180 | } |
| 181 | template <typename Base, typename Exponent> |
| 182 | void int_pow_test() { |
| 183 | Exponent max_exponent = NumTraits<Base>::digits(); |
| 184 | Exponent min_exponent = NumTraits<Exponent>::IsSigned ? -max_exponent : 0; |
| 185 | |
| 186 | for (Exponent exponent = min_exponent; exponent < max_exponent; exponent++) { |
| 187 | test_exponent<Base, Exponent>(exponent); |
| 188 | } |
| 189 | } |
| 190 | |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 191 | template<typename ArrayType> void array(const ArrayType& m) |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 192 | { |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 193 | typedef typename ArrayType::Scalar Scalar; |
Gael Guennebaud | d476cad | 2016-06-25 10:12:06 +0200 | [diff] [blame] | 194 | typedef typename ArrayType::RealScalar RealScalar; |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 195 | typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType; |
| 196 | typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType; |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 197 | |
Hauke Heibel | f1679c7 | 2010-06-20 17:37:56 +0200 | [diff] [blame] | 198 | Index rows = m.rows(); |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 199 | Index cols = m.cols(); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 200 | |
Rasmus Munk Larsen | 98e51c9 | 2022-08-26 16:19:03 +0000 | [diff] [blame] | 201 | ArrayType m1 = ArrayType::Random(rows, cols); |
| 202 | if (NumTraits<RealScalar>::IsInteger && NumTraits<RealScalar>::IsSigned |
| 203 | && !NumTraits<Scalar>::IsComplex) { |
| 204 | // Here we cap the size of the values in m1 such that pow(3)/cube() |
| 205 | // doesn't overflow and result in undefined behavior. Notice that because |
| 206 | // pow(int, int) promotes its inputs and output to double (according to |
| 207 | // the C++ standard), we hvae to make sure that the result fits in 53 bits |
| 208 | // for int64, |
| 209 | RealScalar max_val = |
| 210 | numext::mini(RealScalar(std::cbrt(NumTraits<RealScalar>::highest())), |
| 211 | RealScalar(std::cbrt(1LL << 53)))/2; |
| 212 | m1.array() = (m1.abs().array() <= max_val).select(m1, Scalar(max_val)); |
| 213 | } |
| 214 | ArrayType m2 = ArrayType::Random(rows, cols), |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 215 | m3(rows, cols); |
Christoph Hertzberg | 3be9f5c | 2015-04-16 13:25:20 +0200 | [diff] [blame] | 216 | ArrayType m4 = m1; // copy constructor |
| 217 | VERIFY_IS_APPROX(m1, m4); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 218 | |
Gael Guennebaud | 627595a | 2009-06-10 11:20:30 +0200 | [diff] [blame] | 219 | ColVectorType cv1 = ColVectorType::Random(rows); |
| 220 | RowVectorType rv1 = RowVectorType::Random(cols); |
| 221 | |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 222 | Scalar s1 = internal::random<Scalar>(), |
Hauke Heibel | 8cb3e36 | 2012-03-02 16:27:27 +0100 | [diff] [blame] | 223 | s2 = internal::random<Scalar>(); |
| 224 | |
Gael Guennebaud | 59dc1da | 2008-09-03 17:16:28 +0000 | [diff] [blame] | 225 | // scalar addition |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 226 | VERIFY_IS_APPROX(m1 + s1, s1 + m1); |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 227 | VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1); |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 228 | VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 ); |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 229 | VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1)); |
| 230 | VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1); |
| 231 | VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) ); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 232 | m3 = m1; |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 233 | m3 += s2; |
| 234 | VERIFY_IS_APPROX(m3, m1 + s2); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 235 | m3 = m1; |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 236 | m3 -= s1; |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 237 | VERIFY_IS_APPROX(m3, m1 - s1); |
| 238 | |
Hauke Heibel | f578dc7 | 2010-12-16 17:34:13 +0100 | [diff] [blame] | 239 | // scalar operators via Maps |
Rasmus Munk Larsen | 98e51c9 | 2022-08-26 16:19:03 +0000 | [diff] [blame] | 240 | m3 = m1; m4 = m1; |
| 241 | ArrayType::Map(m4.data(), m4.rows(), m4.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols()); |
| 242 | VERIFY_IS_APPROX(m4, m3 - m2); |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 243 | |
Rasmus Munk Larsen | 98e51c9 | 2022-08-26 16:19:03 +0000 | [diff] [blame] | 244 | m3 = m1; m4 = m1; |
| 245 | ArrayType::Map(m4.data(), m4.rows(), m4.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols()); |
| 246 | VERIFY_IS_APPROX(m4, m3 + m2); |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 247 | |
Rasmus Munk Larsen | 98e51c9 | 2022-08-26 16:19:03 +0000 | [diff] [blame] | 248 | m3 = m1; m4 = m1; |
| 249 | ArrayType::Map(m4.data(), m4.rows(), m4.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols()); |
| 250 | VERIFY_IS_APPROX(m4, m3 * m2); |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 251 | |
Rasmus Munk Larsen | 98e51c9 | 2022-08-26 16:19:03 +0000 | [diff] [blame] | 252 | m3 = m1; m4 = m1; |
Hauke Heibel | f578dc7 | 2010-12-16 17:34:13 +0100 | [diff] [blame] | 253 | m2 = ArrayType::Random(rows,cols); |
| 254 | m2 = (m2==0).select(1,m2); |
Rasmus Munk Larsen | 98e51c9 | 2022-08-26 16:19:03 +0000 | [diff] [blame] | 255 | ArrayType::Map(m4.data(), m4.rows(), m4.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols()); |
| 256 | VERIFY_IS_APPROX(m4, m3 / m2); |
Gael Guennebaud | 05ad083 | 2008-07-19 13:03:23 +0000 | [diff] [blame] | 257 | |
Gael Guennebaud | 59dc1da | 2008-09-03 17:16:28 +0000 | [diff] [blame] | 258 | // reductions |
Gael Guennebaud | 42af587 | 2013-02-23 22:58:14 +0100 | [diff] [blame] | 259 | VERIFY_IS_APPROX(m1.abs().colwise().sum().sum(), m1.abs().sum()); |
| 260 | VERIFY_IS_APPROX(m1.abs().rowwise().sum().sum(), m1.abs().sum()); |
| 261 | using std::abs; |
| 262 | VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.colwise().sum().sum() - m1.sum()), m1.abs().sum()); |
| 263 | VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.rowwise().sum().sum() - m1.sum()), m1.abs().sum()); |
Gael Guennebaud | 28e139a | 2013-02-23 23:06:45 +0100 | [diff] [blame] | 264 | if (!internal::isMuchSmallerThan(abs(m1.sum() - (m1+m2).sum()), m1.abs().sum(), test_precision<Scalar>())) |
Hauke Heibel | 83e3c45 | 2010-12-16 18:53:02 +0100 | [diff] [blame] | 265 | VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum()); |
Gael Guennebaud | 66e99ab | 2016-06-06 15:11:41 +0200 | [diff] [blame] | 266 | VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>())); |
Gael Guennebaud | 627595a | 2009-06-10 11:20:30 +0200 | [diff] [blame] | 267 | |
| 268 | // vector-wise ops |
| 269 | m3 = m1; |
| 270 | VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1); |
| 271 | m3 = m1; |
| 272 | VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1); |
| 273 | m3 = m1; |
| 274 | VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1); |
| 275 | m3 = m1; |
| 276 | VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1); |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 277 | |
Gael Guennebaud | 0a18eec | 2014-09-19 13:25:28 +0200 | [diff] [blame] | 278 | // Conversion from scalar |
| 279 | VERIFY_IS_APPROX((m3 = s1), ArrayType::Constant(rows,cols,s1)); |
| 280 | VERIFY_IS_APPROX((m3 = 1), ArrayType::Constant(rows,cols,1)); |
| 281 | VERIFY_IS_APPROX((m3.topLeftCorner(rows,cols) = 1), ArrayType::Constant(rows,cols,1)); |
| 282 | typedef Array<Scalar, |
| 283 | ArrayType::RowsAtCompileTime==Dynamic?2:ArrayType::RowsAtCompileTime, |
| 284 | ArrayType::ColsAtCompileTime==Dynamic?2:ArrayType::ColsAtCompileTime, |
| 285 | ArrayType::Options> FixedArrayType; |
Gael Guennebaud | 543529d | 2019-01-22 15:30:50 +0100 | [diff] [blame] | 286 | { |
| 287 | FixedArrayType f1(s1); |
| 288 | VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1)); |
| 289 | FixedArrayType f2(numext::real(s1)); |
| 290 | VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1))); |
| 291 | FixedArrayType f3((int)100*numext::real(s1)); |
| 292 | VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1))); |
| 293 | f1.setRandom(); |
| 294 | FixedArrayType f4(f1.data()); |
| 295 | VERIFY_IS_APPROX(f4, f1); |
| 296 | } |
Gael Guennebaud | 543529d | 2019-01-22 15:30:50 +0100 | [diff] [blame] | 297 | { |
| 298 | FixedArrayType f1{s1}; |
| 299 | VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1)); |
| 300 | FixedArrayType f2{numext::real(s1)}; |
| 301 | VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1))); |
| 302 | FixedArrayType f3{(int)100*numext::real(s1)}; |
| 303 | VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1))); |
| 304 | f1.setRandom(); |
| 305 | FixedArrayType f4{f1.data()}; |
| 306 | VERIFY_IS_APPROX(f4, f1); |
| 307 | } |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 308 | |
Gael Guennebaud | d476cad | 2016-06-25 10:12:06 +0200 | [diff] [blame] | 309 | // pow |
| 310 | VERIFY_IS_APPROX(m1.pow(2), m1.square()); |
| 311 | VERIFY_IS_APPROX(pow(m1,2), m1.square()); |
| 312 | VERIFY_IS_APPROX(m1.pow(3), m1.cube()); |
| 313 | VERIFY_IS_APPROX(pow(m1,3), m1.cube()); |
| 314 | VERIFY_IS_APPROX((-m1).pow(3), -m1.cube()); |
| 315 | VERIFY_IS_APPROX(pow(2*m1,3), 8*m1.cube()); |
| 316 | ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2)); |
| 317 | VERIFY_IS_APPROX(Eigen::pow(m1,exponents), m1.square()); |
| 318 | VERIFY_IS_APPROX(m1.pow(exponents), m1.square()); |
| 319 | VERIFY_IS_APPROX(Eigen::pow(2*m1,exponents), 4*m1.square()); |
| 320 | VERIFY_IS_APPROX((2*m1).pow(exponents), 4*m1.square()); |
| 321 | VERIFY_IS_APPROX(Eigen::pow(m1,2*exponents), m1.square().square()); |
| 322 | VERIFY_IS_APPROX(m1.pow(2*exponents), m1.square().square()); |
Gael Guennebaud | e61cee7 | 2016-07-04 11:49:03 +0200 | [diff] [blame] | 323 | VERIFY_IS_APPROX(Eigen::pow(m1(0,0), exponents), ArrayType::Constant(rows,cols,m1(0,0)*m1(0,0))); |
Gael Guennebaud | d476cad | 2016-06-25 10:12:06 +0200 | [diff] [blame] | 324 | |
Gael Guennebaud | 0a18eec | 2014-09-19 13:25:28 +0200 | [diff] [blame] | 325 | // Check possible conflicts with 1D ctor |
| 326 | typedef Array<Scalar, Dynamic, 1> OneDArrayType; |
Gael Guennebaud | 543529d | 2019-01-22 15:30:50 +0100 | [diff] [blame] | 327 | { |
| 328 | OneDArrayType o1(rows); |
| 329 | VERIFY(o1.size()==rows); |
| 330 | OneDArrayType o2(static_cast<int>(rows)); |
| 331 | VERIFY(o2.size()==rows); |
| 332 | } |
Gael Guennebaud | 543529d | 2019-01-22 15:30:50 +0100 | [diff] [blame] | 333 | { |
| 334 | OneDArrayType o1{rows}; |
| 335 | VERIFY(o1.size()==rows); |
| 336 | OneDArrayType o4{int(rows)}; |
| 337 | VERIFY(o4.size()==rows); |
| 338 | } |
Gael Guennebaud | 543529d | 2019-01-22 15:30:50 +0100 | [diff] [blame] | 339 | // Check possible conflicts with 2D ctor |
| 340 | typedef Array<Scalar, Dynamic, Dynamic> TwoDArrayType; |
| 341 | typedef Array<Scalar, 2, 1> ArrayType2; |
| 342 | { |
| 343 | TwoDArrayType o1(rows,cols); |
| 344 | VERIFY(o1.rows()==rows); |
| 345 | VERIFY(o1.cols()==cols); |
| 346 | TwoDArrayType o2(static_cast<int>(rows),static_cast<int>(cols)); |
| 347 | VERIFY(o2.rows()==rows); |
| 348 | VERIFY(o2.cols()==cols); |
| 349 | |
| 350 | ArrayType2 o3(rows,cols); |
| 351 | VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols)); |
| 352 | ArrayType2 o4(static_cast<int>(rows),static_cast<int>(cols)); |
| 353 | VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols)); |
| 354 | } |
Gael Guennebaud | 543529d | 2019-01-22 15:30:50 +0100 | [diff] [blame] | 355 | { |
| 356 | TwoDArrayType o1{rows,cols}; |
| 357 | VERIFY(o1.rows()==rows); |
| 358 | VERIFY(o1.cols()==cols); |
| 359 | TwoDArrayType o2{int(rows),int(cols)}; |
| 360 | VERIFY(o2.rows()==rows); |
| 361 | VERIFY(o2.cols()==cols); |
| 362 | |
| 363 | ArrayType2 o3{rows,cols}; |
| 364 | VERIFY(o3(0)==Scalar(rows) && o3(1)==Scalar(cols)); |
| 365 | ArrayType2 o4{int(rows),int(cols)}; |
| 366 | VERIFY(o4(0)==Scalar(rows) && o4(1)==Scalar(cols)); |
| 367 | } |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 368 | } |
| 369 | |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 370 | template<typename ArrayType> void comparisons(const ArrayType& m) |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 371 | { |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 372 | using std::abs; |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 373 | typedef typename ArrayType::Scalar Scalar; |
Benoit Jacob | 874ff5a | 2009-01-26 17:56:04 +0000 | [diff] [blame] | 374 | typedef typename NumTraits<Scalar>::Real RealScalar; |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 375 | |
Hauke Heibel | f1679c7 | 2010-06-20 17:37:56 +0200 | [diff] [blame] | 376 | Index rows = m.rows(); |
| 377 | Index cols = m.cols(); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 378 | |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 379 | Index r = internal::random<Index>(0, rows-1), |
| 380 | c = internal::random<Index>(0, cols-1); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 381 | |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 382 | ArrayType m1 = ArrayType::Random(rows, cols), |
Gael Guennebaud | 9fc1c92 | 2015-06-22 16:48:27 +0200 | [diff] [blame] | 383 | m2 = ArrayType::Random(rows, cols), |
| 384 | m3(rows, cols), |
| 385 | m4 = m1; |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 386 | |
Gael Guennebaud | 9fc1c92 | 2015-06-22 16:48:27 +0200 | [diff] [blame] | 387 | m4 = (m4.abs()==Scalar(0)).select(1,m4); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 388 | |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 389 | VERIFY(((m1 + Scalar(1)) > m1).all()); |
| 390 | VERIFY(((m1 - Scalar(1)) < m1).all()); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 391 | if (rows*cols>1) |
| 392 | { |
| 393 | m3 = m1; |
| 394 | m3(r,c) += 1; |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 395 | VERIFY(! (m1 < m3).all() ); |
| 396 | VERIFY(! (m1 > m3).all() ); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 397 | } |
Christoph Hertzberg | 494fa99 | 2015-05-07 17:28:40 +0200 | [diff] [blame] | 398 | VERIFY(!(m1 > m2 && m1 < m2).any()); |
| 399 | VERIFY((m1 <= m2 || m1 >= m2).all()); |
Gael Guennebaud | 627595a | 2009-06-10 11:20:30 +0200 | [diff] [blame] | 400 | |
Christoph Hertzberg | 36052c4 | 2013-10-17 15:37:29 +0200 | [diff] [blame] | 401 | // comparisons array to scalar |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 402 | VERIFY( (m1 != (m1(r,c)+1) ).any() ); |
Christoph Hertzberg | 36052c4 | 2013-10-17 15:37:29 +0200 | [diff] [blame] | 403 | VERIFY( (m1 > (m1(r,c)-1) ).any() ); |
| 404 | VERIFY( (m1 < (m1(r,c)+1) ).any() ); |
| 405 | VERIFY( (m1 == m1(r,c) ).any() ); |
| 406 | |
| 407 | // comparisons scalar to array |
| 408 | VERIFY( ( (m1(r,c)+1) != m1).any() ); |
| 409 | VERIFY( ( (m1(r,c)-1) < m1).any() ); |
| 410 | VERIFY( ( (m1(r,c)+1) > m1).any() ); |
| 411 | VERIFY( ( m1(r,c) == m1).any() ); |
Gael Guennebaud | 627595a | 2009-06-10 11:20:30 +0200 | [diff] [blame] | 412 | |
Gael Guennebaud | 59dc1da | 2008-09-03 17:16:28 +0000 | [diff] [blame] | 413 | // test Select |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 414 | VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) ); |
| 415 | VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) ); |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 416 | Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2); |
Gael Guennebaud | 59dc1da | 2008-09-03 17:16:28 +0000 | [diff] [blame] | 417 | for (int j=0; j<cols; ++j) |
| 418 | for (int i=0; i<rows; ++i) |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 419 | m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j); |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 420 | VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid)) |
| 421 | .select(ArrayType::Zero(rows,cols),m1), m3); |
Gael Guennebaud | e14aa8c | 2008-09-03 17:56:06 +0000 | [diff] [blame] | 422 | // shorter versions: |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 423 | VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid)) |
Gael Guennebaud | 59dc1da | 2008-09-03 17:16:28 +0000 | [diff] [blame] | 424 | .select(0,m1), m3); |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 425 | VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid)) |
Gael Guennebaud | 59dc1da | 2008-09-03 17:16:28 +0000 | [diff] [blame] | 426 | .select(m1,0), m3); |
Gael Guennebaud | e14aa8c | 2008-09-03 17:56:06 +0000 | [diff] [blame] | 427 | // even shorter version: |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 428 | VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3); |
Gael Guennebaud | 627595a | 2009-06-10 11:20:30 +0200 | [diff] [blame] | 429 | |
Gael Guennebaud | 56c7e16 | 2009-01-24 15:22:44 +0000 | [diff] [blame] | 430 | // count |
Gael Guennebaud | c70d542 | 2010-01-18 22:54:20 +0100 | [diff] [blame] | 431 | VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols); |
Benoit Jacob | aaaade4 | 2010-05-30 16:00:58 -0400 | [diff] [blame] | 432 | |
Gael Guennebaud | 35c1158 | 2011-05-31 22:17:34 +0200 | [diff] [blame] | 433 | // and/or |
| 434 | VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0); |
| 435 | VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols); |
| 436 | RealScalar a = m1.abs().mean(); |
| 437 | VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count()); |
| 438 | |
Gael Guennebaud | 12e1ebb | 2018-07-12 17:16:40 +0200 | [diff] [blame] | 439 | typedef Array<Index, Dynamic, 1> ArrayOfIndices; |
Gael Guennebaud | 4ebb804 | 2010-06-02 09:45:57 +0200 | [diff] [blame] | 440 | |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 441 | // TODO allows colwise/rowwise for array |
Benoit Jacob | aaaade4 | 2010-05-30 16:00:58 -0400 | [diff] [blame] | 442 | VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose()); |
| 443 | VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols)); |
Benoit Jacob | e800999 | 2008-11-03 22:47:00 +0000 | [diff] [blame] | 444 | } |
| 445 | |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 446 | template<typename ArrayType> void array_real(const ArrayType& m) |
Gael Guennebaud | 0ce5bc0 | 2010-01-27 23:23:59 +0100 | [diff] [blame] | 447 | { |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 448 | using std::abs; |
Gael Guennebaud | 9b9177f | 2013-07-05 13:35:34 +0200 | [diff] [blame] | 449 | using std::sqrt; |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 450 | typedef typename ArrayType::Scalar Scalar; |
Gael Guennebaud | 0ce5bc0 | 2010-01-27 23:23:59 +0100 | [diff] [blame] | 451 | typedef typename NumTraits<Scalar>::Real RealScalar; |
| 452 | |
Hauke Heibel | f1679c7 | 2010-06-20 17:37:56 +0200 | [diff] [blame] | 453 | Index rows = m.rows(); |
| 454 | Index cols = m.cols(); |
Gael Guennebaud | 0ce5bc0 | 2010-01-27 23:23:59 +0100 | [diff] [blame] | 455 | |
Benoit Jacob | e277586 | 2010-04-28 18:51:38 -0400 | [diff] [blame] | 456 | ArrayType m1 = ArrayType::Random(rows, cols), |
Gael Guennebaud | c695cbf | 2013-06-24 13:33:44 +0200 | [diff] [blame] | 457 | m2 = ArrayType::Random(rows, cols), |
Gael Guennebaud | 9fc1c92 | 2015-06-22 16:48:27 +0200 | [diff] [blame] | 458 | m3(rows, cols), |
| 459 | m4 = m1; |
Benoit Steiner | 8314962 | 2015-12-10 13:13:45 -0800 | [diff] [blame] | 460 | |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 461 | m4 = (m4.abs()==Scalar(0)).select(Scalar(1),m4); |
Gael Guennebaud | 0ce5bc0 | 2010-01-27 23:23:59 +0100 | [diff] [blame] | 462 | |
Gael Guennebaud | 9b1ad5e | 2012-03-09 12:08:06 +0100 | [diff] [blame] | 463 | Scalar s1 = internal::random<Scalar>(); |
| 464 | |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 465 | // these tests are mostly to check possible compilation issues with free-functions. |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 466 | VERIFY_IS_APPROX(m1.sin(), sin(m1)); |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 467 | VERIFY_IS_APPROX(m1.cos(), cos(m1)); |
Deanna Hood | f89fcef | 2015-03-11 13:13:30 +1000 | [diff] [blame] | 468 | VERIFY_IS_APPROX(m1.tan(), tan(m1)); |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 469 | VERIFY_IS_APPROX(m1.asin(), asin(m1)); |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 470 | VERIFY_IS_APPROX(m1.acos(), acos(m1)); |
Jitse Niesen | de150b1 | 2014-06-19 15:12:33 +0100 | [diff] [blame] | 471 | VERIFY_IS_APPROX(m1.atan(), atan(m1)); |
Deanna Hood | f89fcef | 2015-03-11 13:13:30 +1000 | [diff] [blame] | 472 | VERIFY_IS_APPROX(m1.sinh(), sinh(m1)); |
| 473 | VERIFY_IS_APPROX(m1.cosh(), cosh(m1)); |
| 474 | VERIFY_IS_APPROX(m1.tanh(), tanh(m1)); |
Rasmus Munk Larsen | 28ba1b2 | 2019-01-11 17:45:37 -0800 | [diff] [blame] | 475 | #if EIGEN_HAS_CXX11_MATH |
| 476 | VERIFY_IS_APPROX(m1.tanh().atanh(), atanh(tanh(m1))); |
| 477 | VERIFY_IS_APPROX(m1.sinh().asinh(), asinh(sinh(m1))); |
| 478 | VERIFY_IS_APPROX(m1.cosh().acosh(), acosh(cosh(m1))); |
| 479 | #endif |
Rasmus Munk Larsen | d6e283b | 2018-08-13 11:14:50 -0700 | [diff] [blame] | 480 | VERIFY_IS_APPROX(m1.logistic(), logistic(m1)); |
Gael Guennebaud | 2f7e261 | 2016-07-08 11:13:55 +0200 | [diff] [blame] | 481 | |
Deanna Hood | a5e4997 | 2015-03-11 08:56:42 +1000 | [diff] [blame] | 482 | VERIFY_IS_APPROX(m1.arg(), arg(m1)); |
Deanna Hood | 31fdd67 | 2015-03-11 06:39:23 +1000 | [diff] [blame] | 483 | VERIFY_IS_APPROX(m1.round(), round(m1)); |
Ilya Tokar | 19876ce | 2019-12-16 16:00:35 -0500 | [diff] [blame] | 484 | VERIFY_IS_APPROX(m1.rint(), rint(m1)); |
Deanna Hood | 31fdd67 | 2015-03-11 06:39:23 +1000 | [diff] [blame] | 485 | VERIFY_IS_APPROX(m1.floor(), floor(m1)); |
| 486 | VERIFY_IS_APPROX(m1.ceil(), ceil(m1)); |
Christoph Hertzberg | 61e0977 | 2015-08-14 17:32:34 +0200 | [diff] [blame] | 487 | VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all()); |
| 488 | VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all()); |
| 489 | VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all()); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 490 | VERIFY_IS_APPROX(m4.inverse(), inverse(m4)); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 491 | VERIFY_IS_APPROX(m1.abs(), abs(m1)); |
| 492 | VERIFY_IS_APPROX(m1.abs2(), abs2(m1)); |
| 493 | VERIFY_IS_APPROX(m1.square(), square(m1)); |
| 494 | VERIFY_IS_APPROX(m1.cube(), cube(m1)); |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 495 | VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval())); |
Gael Guennebaud | 3f32f5e | 2015-11-27 16:27:53 +0100 | [diff] [blame] | 496 | VERIFY_IS_APPROX(m1.sign(), sign(m1)); |
Rasmus Munk Larsen | 6964ae8 | 2020-07-07 01:54:04 +0000 | [diff] [blame] | 497 | VERIFY((m1.sqrt().sign().isNaN() == (Eigen::isnan)(sign(sqrt(m1)))).all()); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 498 | |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 499 | // avoid inf and NaNs so verification doesn't fail |
| 500 | m3 = m4.abs(); |
| 501 | VERIFY_IS_APPROX(m3.sqrt(), sqrt(abs(m3))); |
| 502 | VERIFY_IS_APPROX(m3.rsqrt(), Scalar(1)/sqrt(abs(m3))); |
| 503 | VERIFY_IS_APPROX(rsqrt(m3), Scalar(1)/sqrt(abs(m3))); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 504 | VERIFY_IS_APPROX(m3.log(), log(m3)); |
Gael Guennebaud | 89099b0 | 2016-06-01 17:00:08 +0200 | [diff] [blame] | 505 | VERIFY_IS_APPROX(m3.log1p(), log1p(m3)); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 506 | VERIFY_IS_APPROX(m3.log10(), log10(m3)); |
Rasmus Munk Larsen | f9fac1d | 2020-12-04 21:45:09 +0000 | [diff] [blame] | 507 | VERIFY_IS_APPROX(m3.log2(), log2(m3)); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 508 | |
| 509 | |
| 510 | VERIFY((!(m1>m2) == (m1<=m2)).all()); |
| 511 | |
| 512 | VERIFY_IS_APPROX(sin(m1.asin()), m1); |
| 513 | VERIFY_IS_APPROX(cos(m1.acos()), m1); |
| 514 | VERIFY_IS_APPROX(tan(m1.atan()), m1); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 515 | VERIFY_IS_APPROX(sinh(m1), Scalar(0.5)*(exp(m1)-exp(-m1))); |
| 516 | VERIFY_IS_APPROX(cosh(m1), Scalar(0.5)*(exp(m1)+exp(-m1))); |
| 517 | VERIFY_IS_APPROX(tanh(m1), (Scalar(0.5)*(exp(m1)-exp(-m1)))/(Scalar(0.5)*(exp(m1)+exp(-m1)))); |
| 518 | VERIFY_IS_APPROX(logistic(m1), (Scalar(1)/(Scalar(1)+exp(-m1)))); |
| 519 | VERIFY_IS_APPROX(arg(m1), ((m1<Scalar(0)).template cast<Scalar>())*Scalar(std::acos(Scalar(-1)))); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 520 | VERIFY((round(m1) <= ceil(m1) && round(m1) >= floor(m1)).all()); |
Ilya Tokar | 19876ce | 2019-12-16 16:00:35 -0500 | [diff] [blame] | 521 | VERIFY((rint(m1) <= ceil(m1) && rint(m1) >= floor(m1)).all()); |
| 522 | VERIFY(((ceil(m1) - round(m1)) <= Scalar(0.5) || (round(m1) - floor(m1)) <= Scalar(0.5)).all()); |
| 523 | VERIFY(((ceil(m1) - round(m1)) <= Scalar(1.0) && (round(m1) - floor(m1)) <= Scalar(1.0)).all()); |
| 524 | VERIFY(((ceil(m1) - rint(m1)) <= Scalar(0.5) || (rint(m1) - floor(m1)) <= Scalar(0.5)).all()); |
| 525 | VERIFY(((ceil(m1) - rint(m1)) <= Scalar(1.0) && (rint(m1) - floor(m1)) <= Scalar(1.0)).all()); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 526 | VERIFY((Eigen::isnan)((m1*Scalar(0))/Scalar(0)).all()); |
| 527 | VERIFY((Eigen::isinf)(m4/Scalar(0)).all()); |
| 528 | VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*Scalar(0)/Scalar(0))) && (!(Eigen::isfinite)(m4/Scalar(0)))).all()); |
| 529 | VERIFY_IS_APPROX(inverse(inverse(m4)),m4); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 530 | VERIFY((abs(m1) == m1 || abs(m1) == -m1).all()); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 531 | VERIFY_IS_APPROX(m3, sqrt(abs2(m3))); |
Joel Holdsworth | d5c6657 | 2020-03-19 17:45:20 +0000 | [diff] [blame] | 532 | VERIFY_IS_APPROX(m1.absolute_difference(m2), (m1 > m2).select(m1 - m2, m2 - m1)); |
Gael Guennebaud | 3f32f5e | 2015-11-27 16:27:53 +0100 | [diff] [blame] | 533 | VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() ); |
| 534 | VERIFY_IS_APPROX( m1*m1.sign(),m1.abs()); |
| 535 | VERIFY_IS_APPROX(m1.sign() * m1.abs(), m1); |
Benoit Jacob | 5d63d2c | 2010-04-28 22:42:34 -0400 | [diff] [blame] | 536 | |
Gael Guennebaud | 62670c8 | 2013-06-10 23:40:56 +0200 | [diff] [blame] | 537 | VERIFY_IS_APPROX(numext::abs2(numext::real(m1)) + numext::abs2(numext::imag(m1)), numext::abs2(m1)); |
Gael Guennebaud | 87427d2 | 2019-10-08 09:15:17 +0200 | [diff] [blame] | 538 | VERIFY_IS_APPROX(numext::abs2(Eigen::real(m1)) + numext::abs2(Eigen::imag(m1)), numext::abs2(m1)); |
Benoit Jacob | 5d63d2c | 2010-04-28 22:42:34 -0400 | [diff] [blame] | 539 | if(!NumTraits<Scalar>::IsComplex) |
Gael Guennebaud | 62670c8 | 2013-06-10 23:40:56 +0200 | [diff] [blame] | 540 | VERIFY_IS_APPROX(numext::real(m1), m1); |
Gael Guennebaud | 4ebb804 | 2010-06-02 09:45:57 +0200 | [diff] [blame] | 541 | |
Jitse Niesen | 6fecb6f | 2014-02-24 14:10:17 +0000 | [diff] [blame] | 542 | // shift argument of logarithm so that it is not zero |
| 543 | Scalar smallNumber = NumTraits<Scalar>::dummy_precision(); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 544 | VERIFY_IS_APPROX((m3 + smallNumber).log() , log(abs(m3) + smallNumber)); |
| 545 | VERIFY_IS_APPROX((m3 + smallNumber + Scalar(1)).log() , log1p(abs(m3) + smallNumber)); |
Gael Guennebaud | 4ebb804 | 2010-06-02 09:45:57 +0200 | [diff] [blame] | 546 | |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 547 | VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2)); |
| 548 | VERIFY_IS_APPROX(m1.exp(), exp(m1)); |
| 549 | VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp()); |
Gael Guennebaud | 575ac54 | 2010-06-19 23:17:07 +0200 | [diff] [blame] | 550 | |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 551 | VERIFY_IS_APPROX(m1.expm1(), expm1(m1)); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 552 | VERIFY_IS_APPROX((m3 + smallNumber).exp() - Scalar(1), expm1(abs(m3) + smallNumber)); |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 553 | |
Gael Guennebaud | 575ac54 | 2010-06-19 23:17:07 +0200 | [diff] [blame] | 554 | VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt()); |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 555 | VERIFY_IS_APPROX(pow(m3,RealScalar(0.5)), m3.sqrt()); |
Benoit Steiner | e25e3a0 | 2015-12-03 18:16:35 -0800 | [diff] [blame] | 556 | |
| 557 | VERIFY_IS_APPROX(m3.pow(RealScalar(-0.5)), m3.rsqrt()); |
| 558 | VERIFY_IS_APPROX(pow(m3,RealScalar(-0.5)), m3.rsqrt()); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 559 | |
| 560 | // Avoid inf and NaN. |
| 561 | m3 = (m1.square()<NumTraits<Scalar>::epsilon()).select(Scalar(1),m3); |
| 562 | VERIFY_IS_APPROX(m3.pow(RealScalar(-2)), m3.square().inverse()); |
Rasmus Munk Larsen | cdd8fdc | 2021-01-18 13:25:16 +0000 | [diff] [blame] | 563 | pow_test<Scalar>(); |
Benoit Steiner | e25e3a0 | 2015-12-03 18:16:35 -0800 | [diff] [blame] | 564 | |
Charles Schlosser | e5af9f8 | 2022-08-29 19:23:54 +0000 | [diff] [blame] | 565 | typedef typename internal::make_integer<Scalar>::type SignedInt; |
| 566 | typedef typename std::make_unsigned<SignedInt>::type UnsignedInt; |
| 567 | |
| 568 | int_pow_test<SignedInt, SignedInt>(); |
| 569 | int_pow_test<SignedInt, UnsignedInt>(); |
| 570 | int_pow_test<UnsignedInt, SignedInt>(); |
| 571 | int_pow_test<UnsignedInt, UnsignedInt>(); |
| 572 | |
Antonio Sanchez | 4c42d5e | 2021-01-23 11:54:00 -0800 | [diff] [blame] | 573 | VERIFY_IS_APPROX(log10(m3), log(m3)/numext::log(Scalar(10))); |
| 574 | VERIFY_IS_APPROX(log2(m3), log(m3)/numext::log(Scalar(2))); |
Hauke Heibel | 81c1336 | 2012-03-07 08:58:42 +0100 | [diff] [blame] | 575 | |
| 576 | // scalar by array division |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 577 | const RealScalar tiny = sqrt(std::numeric_limits<RealScalar>::epsilon()); |
Hauke Heibel | dd9365e | 2012-03-09 14:04:13 +0100 | [diff] [blame] | 578 | s1 += Scalar(tiny); |
| 579 | m1 += ArrayType::Constant(rows,cols,Scalar(tiny)); |
Rasmus Munk Larsen | 96dc37a | 2022-01-07 01:10:17 +0000 | [diff] [blame] | 580 | VERIFY_IS_CWISE_APPROX(s1/m1, s1 * m1.inverse()); |
Eugene Brevdo | f736277 | 2015-12-24 21:15:38 -0800 | [diff] [blame] | 581 | |
Gael Guennebaud | c695cbf | 2013-06-24 13:33:44 +0200 | [diff] [blame] | 582 | // check inplace transpose |
| 583 | m3 = m1; |
| 584 | m3.transposeInPlace(); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 585 | VERIFY_IS_APPROX(m3, m1.transpose()); |
Gael Guennebaud | c695cbf | 2013-06-24 13:33:44 +0200 | [diff] [blame] | 586 | m3.transposeInPlace(); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 587 | VERIFY_IS_APPROX(m3, m1); |
Gael Guennebaud | 0ce5bc0 | 2010-01-27 23:23:59 +0100 | [diff] [blame] | 588 | } |
| 589 | |
Jitse Niesen | 837db08 | 2011-05-09 10:17:41 +0100 | [diff] [blame] | 590 | template<typename ArrayType> void array_complex(const ArrayType& m) |
| 591 | { |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 592 | typedef typename ArrayType::Scalar Scalar; |
| 593 | typedef typename NumTraits<Scalar>::Real RealScalar; |
Jitse Niesen | 837db08 | 2011-05-09 10:17:41 +0100 | [diff] [blame] | 594 | |
| 595 | Index rows = m.rows(); |
| 596 | Index cols = m.cols(); |
| 597 | |
| 598 | ArrayType m1 = ArrayType::Random(rows, cols), |
Gael Guennebaud | 9fc1c92 | 2015-06-22 16:48:27 +0200 | [diff] [blame] | 599 | m2(rows, cols), |
| 600 | m4 = m1; |
Srinivas Vasudevan | 218764e | 2016-12-02 14:13:01 -0800 | [diff] [blame] | 601 | |
Gael Guennebaud | 9fc1c92 | 2015-06-22 16:48:27 +0200 | [diff] [blame] | 602 | m4.real() = (m4.real().abs()==RealScalar(0)).select(RealScalar(1),m4.real()); |
| 603 | m4.imag() = (m4.imag().abs()==RealScalar(0)).select(RealScalar(1),m4.imag()); |
Jitse Niesen | 837db08 | 2011-05-09 10:17:41 +0100 | [diff] [blame] | 604 | |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 605 | Array<RealScalar, -1, -1> m3(rows, cols); |
| 606 | |
Jitse Niesen | 837db08 | 2011-05-09 10:17:41 +0100 | [diff] [blame] | 607 | for (Index i = 0; i < m.rows(); ++i) |
| 608 | for (Index j = 0; j < m.cols(); ++j) |
Gael Guennebaud | a76fbbf | 2012-11-06 15:25:50 +0100 | [diff] [blame] | 609 | m2(i,j) = sqrt(m1(i,j)); |
Jitse Niesen | 837db08 | 2011-05-09 10:17:41 +0100 | [diff] [blame] | 610 | |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 611 | // these tests are mostly to check possible compilation issues with free-functions. |
Deanna Hood | f89fcef | 2015-03-11 13:13:30 +1000 | [diff] [blame] | 612 | VERIFY_IS_APPROX(m1.sin(), sin(m1)); |
| 613 | VERIFY_IS_APPROX(m1.cos(), cos(m1)); |
| 614 | VERIFY_IS_APPROX(m1.tan(), tan(m1)); |
| 615 | VERIFY_IS_APPROX(m1.sinh(), sinh(m1)); |
| 616 | VERIFY_IS_APPROX(m1.cosh(), cosh(m1)); |
| 617 | VERIFY_IS_APPROX(m1.tanh(), tanh(m1)); |
Rasmus Munk Larsen | d6e283b | 2018-08-13 11:14:50 -0700 | [diff] [blame] | 618 | VERIFY_IS_APPROX(m1.logistic(), logistic(m1)); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 619 | VERIFY_IS_APPROX(m1.arg(), arg(m1)); |
Christoph Hertzberg | 61e0977 | 2015-08-14 17:32:34 +0200 | [diff] [blame] | 620 | VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all()); |
| 621 | VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all()); |
| 622 | VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all()); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 623 | VERIFY_IS_APPROX(m4.inverse(), inverse(m4)); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 624 | VERIFY_IS_APPROX(m1.log(), log(m1)); |
| 625 | VERIFY_IS_APPROX(m1.log10(), log10(m1)); |
Rasmus Munk Larsen | f9fac1d | 2020-12-04 21:45:09 +0000 | [diff] [blame] | 626 | VERIFY_IS_APPROX(m1.log2(), log2(m1)); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 627 | VERIFY_IS_APPROX(m1.abs(), abs(m1)); |
| 628 | VERIFY_IS_APPROX(m1.abs2(), abs2(m1)); |
| 629 | VERIFY_IS_APPROX(m1.sqrt(), sqrt(m1)); |
| 630 | VERIFY_IS_APPROX(m1.square(), square(m1)); |
| 631 | VERIFY_IS_APPROX(m1.cube(), cube(m1)); |
| 632 | VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval())); |
Gael Guennebaud | 3f32f5e | 2015-11-27 16:27:53 +0100 | [diff] [blame] | 633 | VERIFY_IS_APPROX(m1.sign(), sign(m1)); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 634 | |
| 635 | |
| 636 | VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2)); |
| 637 | VERIFY_IS_APPROX(m1.exp(), exp(m1)); |
| 638 | VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp()); |
| 639 | |
Srinivas Vasudevan | 88062b7 | 2019-12-12 01:56:54 +0000 | [diff] [blame] | 640 | VERIFY_IS_APPROX(m1.expm1(), expm1(m1)); |
| 641 | VERIFY_IS_APPROX(expm1(m1), exp(m1) - 1.); |
| 642 | // Check for larger magnitude complex numbers that expm1 matches exp - 1. |
| 643 | VERIFY_IS_APPROX(expm1(10. * m1), exp(10. * m1) - 1.); |
| 644 | |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 645 | VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1))); |
| 646 | VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1))); |
| 647 | VERIFY_IS_APPROX(tanh(m1), (0.5*(exp(m1)-exp(-m1)))/(0.5*(exp(m1)+exp(-m1)))); |
Rasmus Munk Larsen | d6e283b | 2018-08-13 11:14:50 -0700 | [diff] [blame] | 648 | VERIFY_IS_APPROX(logistic(m1), (1.0/(1.0 + exp(-m1)))); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 649 | |
| 650 | for (Index i = 0; i < m.rows(); ++i) |
| 651 | for (Index j = 0; j < m.cols(); ++j) |
Gael Guennebaud | 87427d2 | 2019-10-08 09:15:17 +0200 | [diff] [blame] | 652 | m3(i,j) = std::atan2(m1(i,j).imag(), m1(i,j).real()); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 653 | VERIFY_IS_APPROX(arg(m1), m3); |
| 654 | |
| 655 | std::complex<RealScalar> zero(0.0,0.0); |
Christoph Hertzberg | 61e0977 | 2015-08-14 17:32:34 +0200 | [diff] [blame] | 656 | VERIFY((Eigen::isnan)(m1*zero/zero).all()); |
Gael Guennebaud | 4a985e7 | 2015-11-20 14:52:08 +0100 | [diff] [blame] | 657 | #if EIGEN_COMP_MSVC |
| 658 | // msvc complex division is not robust |
| 659 | VERIFY((Eigen::isinf)(m4/RealScalar(0)).all()); |
| 660 | #else |
Gael Guennebaud | 40f326e | 2015-06-17 15:33:09 +0200 | [diff] [blame] | 661 | #if EIGEN_COMP_CLANG |
Gael Guennebaud | 4a985e7 | 2015-11-20 14:52:08 +0100 | [diff] [blame] | 662 | // clang's complex division is notoriously broken too |
Christoph Hertzberg | 61e0977 | 2015-08-14 17:32:34 +0200 | [diff] [blame] | 663 | if((numext::isinf)(m4(0,0)/RealScalar(0))) { |
Gael Guennebaud | 40f326e | 2015-06-17 15:33:09 +0200 | [diff] [blame] | 664 | #endif |
Gael Guennebaud | 4a985e7 | 2015-11-20 14:52:08 +0100 | [diff] [blame] | 665 | VERIFY((Eigen::isinf)(m4/zero).all()); |
Gael Guennebaud | 40f326e | 2015-06-17 15:33:09 +0200 | [diff] [blame] | 666 | #if EIGEN_COMP_CLANG |
| 667 | } |
| 668 | else |
| 669 | { |
Christoph Hertzberg | 61e0977 | 2015-08-14 17:32:34 +0200 | [diff] [blame] | 670 | VERIFY((Eigen::isinf)(m4.real()/zero.real()).all()); |
Gael Guennebaud | 40f326e | 2015-06-17 15:33:09 +0200 | [diff] [blame] | 671 | } |
| 672 | #endif |
Gael Guennebaud | 4a985e7 | 2015-11-20 14:52:08 +0100 | [diff] [blame] | 673 | #endif // MSVC |
| 674 | |
Christoph Hertzberg | 61e0977 | 2015-08-14 17:32:34 +0200 | [diff] [blame] | 675 | VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*zero/zero)) && (!(Eigen::isfinite)(m1/zero))).all()); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 676 | |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 677 | VERIFY_IS_APPROX(inverse(inverse(m4)),m4); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 678 | VERIFY_IS_APPROX(conj(m1.conjugate()), m1); |
Gael Guennebaud | 87427d2 | 2019-10-08 09:15:17 +0200 | [diff] [blame] | 679 | VERIFY_IS_APPROX(abs(m1), sqrt(square(m1.real())+square(m1.imag()))); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 680 | VERIFY_IS_APPROX(abs(m1), sqrt(abs2(m1))); |
| 681 | VERIFY_IS_APPROX(log10(m1), log(m1)/log(10)); |
Rasmus Munk Larsen | f9fac1d | 2020-12-04 21:45:09 +0000 | [diff] [blame] | 682 | VERIFY_IS_APPROX(log2(m1), log(m1)/log(2)); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 683 | |
Gael Guennebaud | 3f32f5e | 2015-11-27 16:27:53 +0100 | [diff] [blame] | 684 | VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() ); |
| 685 | VERIFY_IS_APPROX( m1.sign() * m1.abs(), m1); |
| 686 | |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 687 | // scalar by array division |
Eugene Brevdo | f736277 | 2015-12-24 21:15:38 -0800 | [diff] [blame] | 688 | Scalar s1 = internal::random<Scalar>(); |
Benoit Steiner | 2d74ef9 | 2016-05-17 07:20:11 -0700 | [diff] [blame] | 689 | const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon()); |
Deanna Hood | 41b717d | 2015-03-18 03:11:03 +1000 | [diff] [blame] | 690 | s1 += Scalar(tiny); |
| 691 | m1 += ArrayType::Constant(rows,cols,Scalar(tiny)); |
| 692 | VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse()); |
| 693 | |
| 694 | // check inplace transpose |
| 695 | m2 = m1; |
| 696 | m2.transposeInPlace(); |
| 697 | VERIFY_IS_APPROX(m2, m1.transpose()); |
| 698 | m2.transposeInPlace(); |
| 699 | VERIFY_IS_APPROX(m2, m1); |
Rasmus Munk Larsen | b47c777 | 2020-04-27 18:55:15 -0700 | [diff] [blame] | 700 | // Check vectorized inplace transpose. |
Rasmus Munk Larsen | 74ec8e6 | 2020-05-07 17:29:56 +0000 | [diff] [blame] | 701 | ArrayType m5 = ArrayType::Random(131, 131); |
Rasmus Munk Larsen | b47c777 | 2020-04-27 18:55:15 -0700 | [diff] [blame] | 702 | ArrayType m6 = m5; |
| 703 | m6.transposeInPlace(); |
| 704 | VERIFY_IS_APPROX(m6, m5.transpose()); |
Jitse Niesen | 837db08 | 2011-05-09 10:17:41 +0100 | [diff] [blame] | 705 | } |
| 706 | |
Abraham Bachrach | 039408c | 2012-01-11 11:00:30 -0500 | [diff] [blame] | 707 | template<typename ArrayType> void min_max(const ArrayType& m) |
| 708 | { |
Abraham Bachrach | 039408c | 2012-01-11 11:00:30 -0500 | [diff] [blame] | 709 | typedef typename ArrayType::Scalar Scalar; |
| 710 | |
| 711 | Index rows = m.rows(); |
| 712 | Index cols = m.cols(); |
| 713 | |
| 714 | ArrayType m1 = ArrayType::Random(rows, cols); |
| 715 | |
| 716 | // min/max with array |
| 717 | Scalar maxM1 = m1.maxCoeff(); |
| 718 | Scalar minM1 = m1.minCoeff(); |
| 719 | |
| 720 | VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1))); |
| 721 | VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1))); |
| 722 | |
| 723 | VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1))); |
| 724 | VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1))); |
| 725 | |
| 726 | // min/max with scalar input |
| 727 | VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1)); |
| 728 | VERIFY_IS_APPROX(m1, (m1.min)( maxM1)); |
| 729 | |
| 730 | VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1)); |
| 731 | VERIFY_IS_APPROX(m1, (m1.max)( minM1)); |
| 732 | |
Rasmus Munk Larsen | 841c898 | 2021-02-24 17:49:20 -0800 | [diff] [blame] | 733 | |
| 734 | // min/max with various NaN propagation options. |
| 735 | if (m1.size() > 1 && !NumTraits<Scalar>::IsInteger) { |
Antonio Sanchez | 8dfe102 | 2021-03-16 20:12:46 -0700 | [diff] [blame] | 736 | m1(0,0) = NumTraits<Scalar>::quiet_NaN(); |
Rasmus Munk Larsen | 841c898 | 2021-02-24 17:49:20 -0800 | [diff] [blame] | 737 | maxM1 = m1.template maxCoeff<PropagateNaN>(); |
| 738 | minM1 = m1.template minCoeff<PropagateNaN>(); |
| 739 | VERIFY((numext::isnan)(maxM1)); |
| 740 | VERIFY((numext::isnan)(minM1)); |
| 741 | |
| 742 | maxM1 = m1.template maxCoeff<PropagateNumbers>(); |
| 743 | minM1 = m1.template minCoeff<PropagateNumbers>(); |
| 744 | VERIFY(!(numext::isnan)(maxM1)); |
| 745 | VERIFY(!(numext::isnan)(minM1)); |
| 746 | } |
Abraham Bachrach | 039408c | 2012-01-11 11:00:30 -0500 | [diff] [blame] | 747 | } |
| 748 | |
Antonio Sanchez | fc9d352 | 2021-08-17 20:04:48 -0700 | [diff] [blame] | 749 | template<int N> |
| 750 | struct shift_left { |
| 751 | template<typename Scalar> |
| 752 | Scalar operator()(const Scalar& v) const { |
| 753 | return v << N; |
| 754 | } |
| 755 | }; |
| 756 | |
| 757 | template<int N> |
| 758 | struct arithmetic_shift_right { |
| 759 | template<typename Scalar> |
| 760 | Scalar operator()(const Scalar& v) const { |
| 761 | return v >> N; |
| 762 | } |
| 763 | }; |
| 764 | |
| 765 | template<typename ArrayType> void array_integer(const ArrayType& m) |
| 766 | { |
| 767 | Index rows = m.rows(); |
| 768 | Index cols = m.cols(); |
| 769 | |
| 770 | ArrayType m1 = ArrayType::Random(rows, cols), |
| 771 | m2(rows, cols); |
| 772 | |
| 773 | m2 = m1.template shiftLeft<2>(); |
| 774 | VERIFY( (m2 == m1.unaryExpr(shift_left<2>())).all() ); |
| 775 | m2 = m1.template shiftLeft<9>(); |
| 776 | VERIFY( (m2 == m1.unaryExpr(shift_left<9>())).all() ); |
| 777 | |
| 778 | m2 = m1.template shiftRight<2>(); |
| 779 | VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<2>())).all() ); |
| 780 | m2 = m1.template shiftRight<9>(); |
| 781 | VERIFY( (m2 == m1.unaryExpr(arithmetic_shift_right<9>())).all() ); |
| 782 | } |
| 783 | |
Gael Guennebaud | e0f6d35 | 2018-09-20 18:07:32 +0200 | [diff] [blame] | 784 | EIGEN_DECLARE_TEST(array_cwise) |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 785 | { |
| 786 | for(int i = 0; i < g_repeat; i++) { |
Hauke Heibel | 2d0dfe5 | 2010-12-16 17:36:10 +0100 | [diff] [blame] | 787 | CALL_SUBTEST_1( array(Array<float, 1, 1>()) ); |
Hauke Heibel | b31e124 | 2010-12-16 19:07:23 +0100 | [diff] [blame] | 788 | CALL_SUBTEST_2( array(Array22f()) ); |
| 789 | CALL_SUBTEST_3( array(Array44d()) ); |
Gael Guennebaud | a8f66fe | 2011-07-12 14:41:00 +0200 | [diff] [blame] | 790 | CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
| 791 | CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
| 792 | CALL_SUBTEST_6( array(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
Gael Guennebaud | 543529d | 2019-01-22 15:30:50 +0100 | [diff] [blame] | 793 | 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 Sanchez | fc9d352 | 2021-08-17 20:04:48 -0700 | [diff] [blame] | 794 | CALL_SUBTEST_6( array_integer(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
| 795 | 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 Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 796 | } |
Hauke Heibel | 2d0dfe5 | 2010-12-16 17:36:10 +0100 | [diff] [blame] | 797 | for(int i = 0; i < g_repeat; i++) { |
| 798 | CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) ); |
Hauke Heibel | b31e124 | 2010-12-16 19:07:23 +0100 | [diff] [blame] | 799 | CALL_SUBTEST_2( comparisons(Array22f()) ); |
| 800 | CALL_SUBTEST_3( comparisons(Array44d()) ); |
Gael Guennebaud | a8f66fe | 2011-07-12 14:41:00 +0200 | [diff] [blame] | 801 | CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
| 802 | CALL_SUBTEST_6( comparisons(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
Hauke Heibel | 2d0dfe5 | 2010-12-16 17:36:10 +0100 | [diff] [blame] | 803 | } |
| 804 | for(int i = 0; i < g_repeat; i++) { |
Abraham Bachrach | 039408c | 2012-01-11 11:00:30 -0500 | [diff] [blame] | 805 | CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) ); |
| 806 | CALL_SUBTEST_2( min_max(Array22f()) ); |
| 807 | CALL_SUBTEST_3( min_max(Array44d()) ); |
| 808 | CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
| 809 | CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
| 810 | } |
| 811 | for(int i = 0; i < g_repeat; i++) { |
Hauke Heibel | 2d0dfe5 | 2010-12-16 17:36:10 +0100 | [diff] [blame] | 812 | CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) ); |
Hauke Heibel | b31e124 | 2010-12-16 19:07:23 +0100 | [diff] [blame] | 813 | CALL_SUBTEST_2( array_real(Array22f()) ); |
| 814 | CALL_SUBTEST_3( array_real(Array44d()) ); |
Gael Guennebaud | a8f66fe | 2011-07-12 14:41:00 +0200 | [diff] [blame] | 815 | CALL_SUBTEST_5( array_real(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
Antonio Sanchez | f0e46ed | 2021-01-22 11:10:54 -0800 | [diff] [blame] | 816 | CALL_SUBTEST_7( array_real(Array<Eigen::half, 32, 32>()) ); |
| 817 | CALL_SUBTEST_8( array_real(Array<Eigen::bfloat16, 32, 32>()) ); |
Hauke Heibel | 2d0dfe5 | 2010-12-16 17:36:10 +0100 | [diff] [blame] | 818 | } |
Jitse Niesen | 837db08 | 2011-05-09 10:17:41 +0100 | [diff] [blame] | 819 | for(int i = 0; i < g_repeat; i++) { |
Gael Guennebaud | a8f66fe | 2011-07-12 14:41:00 +0200 | [diff] [blame] | 820 | CALL_SUBTEST_4( array_complex(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); |
Jitse Niesen | 837db08 | 2011-05-09 10:17:41 +0100 | [diff] [blame] | 821 | } |
Benoit Jacob | 5d63d2c | 2010-04-28 22:42:34 -0400 | [diff] [blame] | 822 | |
Hauke Heibel | 2d0dfe5 | 2010-12-16 17:36:10 +0100 | [diff] [blame] | 823 | VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value)); |
| 824 | VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value)); |
| 825 | VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value)); |
Gael Guennebaud | 64fcfd3 | 2016-06-14 11:26:57 +0200 | [diff] [blame] | 826 | typedef CwiseUnaryOp<internal::scalar_abs_op<double>, ArrayXd > Xpr; |
Hauke Heibel | 2d0dfe5 | 2010-12-16 17:36:10 +0100 | [diff] [blame] | 827 | VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type, |
| 828 | ArrayBase<Xpr> |
| 829 | >::value)); |
Gael Guennebaud | 8cef541 | 2008-06-21 17:28:07 +0000 | [diff] [blame] | 830 | } |