1 // This file is part of Eigen, a lightweight C++ template library 2 // for linear algebra. 3 // 4 // Copyright (C) 2008-2009 Gael Guennebaud <gael.guennebaud (at) inria.fr> 5 // 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/. 9 10 #include "main.h" 11 12 template<typename ArrayType> void array(const ArrayType& m) 13 { 14 typedef typename ArrayType::Index Index; 15 typedef typename ArrayType::Scalar Scalar; 16 typedef typename NumTraits<Scalar>::Real RealScalar; 17 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType; 18 typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType; 19 20 Index rows = m.rows(); 21 Index cols = m.cols(); 22 23 ArrayType m1 = ArrayType::Random(rows, cols), 24 m2 = ArrayType::Random(rows, cols), 25 m3(rows, cols); 26 27 ColVectorType cv1 = ColVectorType::Random(rows); 28 RowVectorType rv1 = RowVectorType::Random(cols); 29 30 Scalar s1 = internal::random<Scalar>(), 31 s2 = internal::random<Scalar>(); 32 33 // scalar addition 34 VERIFY_IS_APPROX(m1 + s1, s1 + m1); 35 VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1); 36 VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 ); 37 VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1)); 38 VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1); 39 VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) ); 40 m3 = m1; 41 m3 += s2; 42 VERIFY_IS_APPROX(m3, m1 + s2); 43 m3 = m1; 44 m3 -= s1; 45 VERIFY_IS_APPROX(m3, m1 - s1); 46 47 // scalar operators via Maps 48 m3 = m1; 49 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols()); 50 VERIFY_IS_APPROX(m1, m3 - m2); 51 52 m3 = m1; 53 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols()); 54 VERIFY_IS_APPROX(m1, m3 + m2); 55 56 m3 = m1; 57 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols()); 58 VERIFY_IS_APPROX(m1, m3 * m2); 59 60 m3 = m1; 61 m2 = ArrayType::Random(rows,cols); 62 m2 = (m2==0).select(1,m2); 63 ArrayType::Map(m1.data(), m1.rows(), m1.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols()); 64 VERIFY_IS_APPROX(m1, m3 / m2); 65 66 // reductions 67 VERIFY_IS_APPROX(m1.colwise().sum().sum(), m1.sum()); 68 VERIFY_IS_APPROX(m1.rowwise().sum().sum(), m1.sum()); 69 if (!internal::isApprox(m1.sum(), (m1+m2).sum(), test_precision<Scalar>())) 70 VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum()); 71 VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar>())); 72 73 // vector-wise ops 74 m3 = m1; 75 VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1); 76 m3 = m1; 77 VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1); 78 m3 = m1; 79 VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1); 80 m3 = m1; 81 VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1); 82 } 83 84 template<typename ArrayType> void comparisons(const ArrayType& m) 85 { 86 typedef typename ArrayType::Index Index; 87 typedef typename ArrayType::Scalar Scalar; 88 typedef typename NumTraits<Scalar>::Real RealScalar; 89 typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> VectorType; 90 91 Index rows = m.rows(); 92 Index cols = m.cols(); 93 94 Index r = internal::random<Index>(0, rows-1), 95 c = internal::random<Index>(0, cols-1); 96 97 ArrayType m1 = ArrayType::Random(rows, cols), 98 m2 = ArrayType::Random(rows, cols), 99 m3(rows, cols); 100 101 VERIFY(((m1 + Scalar(1)) > m1).all()); 102 VERIFY(((m1 - Scalar(1)) < m1).all()); 103 if (rows*cols>1) 104 { 105 m3 = m1; 106 m3(r,c) += 1; 107 VERIFY(! (m1 < m3).all() ); 108 VERIFY(! (m1 > m3).all() ); 109 } 110 111 // comparisons to scalar 112 VERIFY( (m1 != (m1(r,c)+1) ).any() ); 113 VERIFY( (m1 > (m1(r,c)-1) ).any() ); 114 VERIFY( (m1 < (m1(r,c)+1) ).any() ); 115 VERIFY( (m1 == m1(r,c) ).any() ); 116 117 // test Select 118 VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) ); 119 VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) ); 120 Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2); 121 for (int j=0; j<cols; ++j) 122 for (int i=0; i<rows; ++i) 123 m3(i,j) = internal::abs(m1(i,j))<mid ? 0 : m1(i,j); 124 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid)) 125 .select(ArrayType::Zero(rows,cols),m1), m3); 126 // shorter versions: 127 VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid)) 128 .select(0,m1), m3); 129 VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid)) 130 .select(m1,0), m3); 131 // even shorter version: 132 VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3); 133 134 // count 135 VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols); 136 137 // and/or 138 VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0); 139 VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols); 140 RealScalar a = m1.abs().mean(); 141 VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count()); 142 143 typedef Array<typename ArrayType::Index, Dynamic, 1> ArrayOfIndices; 144 145 // TODO allows colwise/rowwise for array 146 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose()); 147 VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols)); 148 } 149 150 template<typename ArrayType> void array_real(const ArrayType& m) 151 { 152 typedef typename ArrayType::Index Index; 153 typedef typename ArrayType::Scalar Scalar; 154 typedef typename NumTraits<Scalar>::Real RealScalar; 155 156 Index rows = m.rows(); 157 Index cols = m.cols(); 158 159 ArrayType m1 = ArrayType::Random(rows, cols), 160 m2 = ArrayType::Random(rows, cols), 161 m3(rows, cols); 162 163 Scalar s1 = internal::random<Scalar>(); 164 165 // these tests are mostly to check possible compilation issues. 166 VERIFY_IS_APPROX(m1.sin(), std::sin(m1)); 167 VERIFY_IS_APPROX(m1.sin(), internal::sin(m1)); 168 VERIFY_IS_APPROX(m1.cos(), std::cos(m1)); 169 VERIFY_IS_APPROX(m1.cos(), internal::cos(m1)); 170 VERIFY_IS_APPROX(m1.asin(), std::asin(m1)); 171 VERIFY_IS_APPROX(m1.asin(), internal::asin(m1)); 172 VERIFY_IS_APPROX(m1.acos(), std::acos(m1)); 173 VERIFY_IS_APPROX(m1.acos(), internal::acos(m1)); 174 VERIFY_IS_APPROX(m1.tan(), std::tan(m1)); 175 VERIFY_IS_APPROX(m1.tan(), internal::tan(m1)); 176 177 VERIFY_IS_APPROX(internal::cos(m1+RealScalar(3)*m2), internal::cos((m1+RealScalar(3)*m2).eval())); 178 VERIFY_IS_APPROX(std::cos(m1+RealScalar(3)*m2), std::cos((m1+RealScalar(3)*m2).eval())); 179 180 VERIFY_IS_APPROX(m1.abs().sqrt(), std::sqrt(std::abs(m1))); 181 VERIFY_IS_APPROX(m1.abs().sqrt(), internal::sqrt(internal::abs(m1))); 182 VERIFY_IS_APPROX(m1.abs(), internal::sqrt(internal::abs2(m1))); 183 184 VERIFY_IS_APPROX(internal::abs2(internal::real(m1)) + internal::abs2(internal::imag(m1)), internal::abs2(m1)); 185 VERIFY_IS_APPROX(internal::abs2(std::real(m1)) + internal::abs2(std::imag(m1)), internal::abs2(m1)); 186 if(!NumTraits<Scalar>::IsComplex) 187 VERIFY_IS_APPROX(internal::real(m1), m1); 188 189 VERIFY_IS_APPROX(m1.abs().log(), std::log(std::abs(m1))); 190 VERIFY_IS_APPROX(m1.abs().log(), internal::log(internal::abs(m1))); 191 192 VERIFY_IS_APPROX(m1.exp(), std::exp(m1)); 193 VERIFY_IS_APPROX(m1.exp() * m2.exp(), std::exp(m1+m2)); 194 VERIFY_IS_APPROX(m1.exp(), internal::exp(m1)); 195 VERIFY_IS_APPROX(m1.exp() / m2.exp(), std::exp(m1-m2)); 196 197 VERIFY_IS_APPROX(m1.pow(2), m1.square()); 198 VERIFY_IS_APPROX(std::pow(m1,2), m1.square()); 199 200 ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2)); 201 VERIFY_IS_APPROX(std::pow(m1,exponents), m1.square()); 202 203 m3 = m1.abs(); 204 VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt()); 205 VERIFY_IS_APPROX(std::pow(m3,RealScalar(0.5)), m3.sqrt()); 206 207 // scalar by array division 208 const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon()); 209 s1 += Scalar(tiny); 210 m1 += ArrayType::Constant(rows,cols,Scalar(tiny)); 211 VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse()); 212 } 213 214 template<typename ArrayType> void array_complex(const ArrayType& m) 215 { 216 typedef typename ArrayType::Index Index; 217 218 Index rows = m.rows(); 219 Index cols = m.cols(); 220 221 ArrayType m1 = ArrayType::Random(rows, cols), 222 m2(rows, cols); 223 224 for (Index i = 0; i < m.rows(); ++i) 225 for (Index j = 0; j < m.cols(); ++j) 226 m2(i,j) = std::sqrt(m1(i,j)); 227 228 VERIFY_IS_APPROX(m1.sqrt(), m2); 229 VERIFY_IS_APPROX(m1.sqrt(), std::sqrt(m1)); 230 VERIFY_IS_APPROX(m1.sqrt(), internal::sqrt(m1)); 231 } 232 233 template<typename ArrayType> void min_max(const ArrayType& m) 234 { 235 typedef typename ArrayType::Index Index; 236 typedef typename ArrayType::Scalar Scalar; 237 238 Index rows = m.rows(); 239 Index cols = m.cols(); 240 241 ArrayType m1 = ArrayType::Random(rows, cols); 242 243 // min/max with array 244 Scalar maxM1 = m1.maxCoeff(); 245 Scalar minM1 = m1.minCoeff(); 246 247 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1))); 248 VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1))); 249 250 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1))); 251 VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1))); 252 253 // min/max with scalar input 254 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1)); 255 VERIFY_IS_APPROX(m1, (m1.min)( maxM1)); 256 257 VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1)); 258 VERIFY_IS_APPROX(m1, (m1.max)( minM1)); 259 260 } 261 262 void test_array() 263 { 264 for(int i = 0; i < g_repeat; i++) { 265 CALL_SUBTEST_1( array(Array<float, 1, 1>()) ); 266 CALL_SUBTEST_2( array(Array22f()) ); 267 CALL_SUBTEST_3( array(Array44d()) ); 268 CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 269 CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 270 CALL_SUBTEST_6( array(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 271 } 272 for(int i = 0; i < g_repeat; i++) { 273 CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) ); 274 CALL_SUBTEST_2( comparisons(Array22f()) ); 275 CALL_SUBTEST_3( comparisons(Array44d()) ); 276 CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 277 CALL_SUBTEST_6( comparisons(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 278 } 279 for(int i = 0; i < g_repeat; i++) { 280 CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) ); 281 CALL_SUBTEST_2( min_max(Array22f()) ); 282 CALL_SUBTEST_3( min_max(Array44d()) ); 283 CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 284 CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 285 } 286 for(int i = 0; i < g_repeat; i++) { 287 CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) ); 288 CALL_SUBTEST_2( array_real(Array22f()) ); 289 CALL_SUBTEST_3( array_real(Array44d()) ); 290 CALL_SUBTEST_5( array_real(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 291 } 292 for(int i = 0; i < g_repeat; i++) { 293 CALL_SUBTEST_4( array_complex(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) ); 294 } 295 296 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value)); 297 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value)); 298 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value)); 299 typedef CwiseUnaryOp<internal::scalar_sum_op<double>, ArrayXd > Xpr; 300 VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type, 301 ArrayBase<Xpr> 302 >::value)); 303 } 304