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