• Home
  • History
  • Annotate
  • 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 typename ArrayType::RealScalar 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   ArrayType m4 = m1; // copy constructor
27   VERIFY_IS_APPROX(m1, m4);
28 
29   ColVectorType cv1 = ColVectorType::Random(rows);
30   RowVectorType rv1 = RowVectorType::Random(cols);
31 
32   Scalar  s1 = internal::random<Scalar>(),
33           s2 = internal::random<Scalar>();
34 
35   // scalar addition
36   VERIFY_IS_APPROX(m1 + s1, s1 + m1);
37   VERIFY_IS_APPROX(m1 + s1, ArrayType::Constant(rows,cols,s1) + m1);
38   VERIFY_IS_APPROX(s1 - m1, (-m1)+s1 );
39   VERIFY_IS_APPROX(m1 - s1, m1 - ArrayType::Constant(rows,cols,s1));
40   VERIFY_IS_APPROX(s1 - m1, ArrayType::Constant(rows,cols,s1) - m1);
41   VERIFY_IS_APPROX((m1*Scalar(2)) - s2, (m1+m1) - ArrayType::Constant(rows,cols,s2) );
42   m3 = m1;
43   m3 += s2;
44   VERIFY_IS_APPROX(m3, m1 + s2);
45   m3 = m1;
46   m3 -= s1;
47   VERIFY_IS_APPROX(m3, m1 - s1);
48 
49   // scalar operators via Maps
50   m3 = m1;
51   ArrayType::Map(m1.data(), m1.rows(), m1.cols()) -= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
52   VERIFY_IS_APPROX(m1, m3 - m2);
53 
54   m3 = m1;
55   ArrayType::Map(m1.data(), m1.rows(), m1.cols()) += ArrayType::Map(m2.data(), m2.rows(), m2.cols());
56   VERIFY_IS_APPROX(m1, m3 + m2);
57 
58   m3 = m1;
59   ArrayType::Map(m1.data(), m1.rows(), m1.cols()) *= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
60   VERIFY_IS_APPROX(m1, m3 * m2);
61 
62   m3 = m1;
63   m2 = ArrayType::Random(rows,cols);
64   m2 = (m2==0).select(1,m2);
65   ArrayType::Map(m1.data(), m1.rows(), m1.cols()) /= ArrayType::Map(m2.data(), m2.rows(), m2.cols());
66   VERIFY_IS_APPROX(m1, m3 / m2);
67 
68   // reductions
69   VERIFY_IS_APPROX(m1.abs().colwise().sum().sum(), m1.abs().sum());
70   VERIFY_IS_APPROX(m1.abs().rowwise().sum().sum(), m1.abs().sum());
71   using std::abs;
72   VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.colwise().sum().sum() - m1.sum()), m1.abs().sum());
73   VERIFY_IS_MUCH_SMALLER_THAN(abs(m1.rowwise().sum().sum() - m1.sum()), m1.abs().sum());
74   if (!internal::isMuchSmallerThan(abs(m1.sum() - (m1+m2).sum()), m1.abs().sum(), test_precision<Scalar>()))
75       VERIFY_IS_NOT_APPROX(((m1+m2).rowwise().sum()).sum(), m1.sum());
76   VERIFY_IS_APPROX(m1.colwise().sum(), m1.colwise().redux(internal::scalar_sum_op<Scalar,Scalar>()));
77 
78   // vector-wise ops
79   m3 = m1;
80   VERIFY_IS_APPROX(m3.colwise() += cv1, m1.colwise() + cv1);
81   m3 = m1;
82   VERIFY_IS_APPROX(m3.colwise() -= cv1, m1.colwise() - cv1);
83   m3 = m1;
84   VERIFY_IS_APPROX(m3.rowwise() += rv1, m1.rowwise() + rv1);
85   m3 = m1;
86   VERIFY_IS_APPROX(m3.rowwise() -= rv1, m1.rowwise() - rv1);
87 
88   // Conversion from scalar
89   VERIFY_IS_APPROX((m3 = s1), ArrayType::Constant(rows,cols,s1));
90   VERIFY_IS_APPROX((m3 = 1),  ArrayType::Constant(rows,cols,1));
91   VERIFY_IS_APPROX((m3.topLeftCorner(rows,cols) = 1),  ArrayType::Constant(rows,cols,1));
92   typedef Array<Scalar,
93                 ArrayType::RowsAtCompileTime==Dynamic?2:ArrayType::RowsAtCompileTime,
94                 ArrayType::ColsAtCompileTime==Dynamic?2:ArrayType::ColsAtCompileTime,
95                 ArrayType::Options> FixedArrayType;
96   FixedArrayType f1(s1);
97   VERIFY_IS_APPROX(f1, FixedArrayType::Constant(s1));
98   FixedArrayType f2(numext::real(s1));
99   VERIFY_IS_APPROX(f2, FixedArrayType::Constant(numext::real(s1)));
100   FixedArrayType f3((int)100*numext::real(s1));
101   VERIFY_IS_APPROX(f3, FixedArrayType::Constant((int)100*numext::real(s1)));
102   f1.setRandom();
103   FixedArrayType f4(f1.data());
104   VERIFY_IS_APPROX(f4, f1);
105 
106   // pow
107   VERIFY_IS_APPROX(m1.pow(2), m1.square());
108   VERIFY_IS_APPROX(pow(m1,2), m1.square());
109   VERIFY_IS_APPROX(m1.pow(3), m1.cube());
110   VERIFY_IS_APPROX(pow(m1,3), m1.cube());
111   VERIFY_IS_APPROX((-m1).pow(3), -m1.cube());
112   VERIFY_IS_APPROX(pow(2*m1,3), 8*m1.cube());
113   ArrayType exponents = ArrayType::Constant(rows, cols, RealScalar(2));
114   VERIFY_IS_APPROX(Eigen::pow(m1,exponents), m1.square());
115   VERIFY_IS_APPROX(m1.pow(exponents), m1.square());
116   VERIFY_IS_APPROX(Eigen::pow(2*m1,exponents), 4*m1.square());
117   VERIFY_IS_APPROX((2*m1).pow(exponents), 4*m1.square());
118   VERIFY_IS_APPROX(Eigen::pow(m1,2*exponents), m1.square().square());
119   VERIFY_IS_APPROX(m1.pow(2*exponents), m1.square().square());
120   VERIFY_IS_APPROX(Eigen::pow(m1(0,0), exponents), ArrayType::Constant(rows,cols,m1(0,0)*m1(0,0)));
121 
122   // Check possible conflicts with 1D ctor
123   typedef Array<Scalar, Dynamic, 1> OneDArrayType;
124   OneDArrayType o1(rows);
125   VERIFY(o1.size()==rows);
126   OneDArrayType o4((int)rows);
127   VERIFY(o4.size()==rows);
128 }
129 
comparisons(const ArrayType & m)130 template<typename ArrayType> void comparisons(const ArrayType& m)
131 {
132   using std::abs;
133   typedef typename ArrayType::Index Index;
134   typedef typename ArrayType::Scalar Scalar;
135   typedef typename NumTraits<Scalar>::Real RealScalar;
136 
137   Index rows = m.rows();
138   Index cols = m.cols();
139 
140   Index r = internal::random<Index>(0, rows-1),
141         c = internal::random<Index>(0, cols-1);
142 
143   ArrayType m1 = ArrayType::Random(rows, cols),
144             m2 = ArrayType::Random(rows, cols),
145             m3(rows, cols),
146             m4 = m1;
147 
148   m4 = (m4.abs()==Scalar(0)).select(1,m4);
149 
150   VERIFY(((m1 + Scalar(1)) > m1).all());
151   VERIFY(((m1 - Scalar(1)) < m1).all());
152   if (rows*cols>1)
153   {
154     m3 = m1;
155     m3(r,c) += 1;
156     VERIFY(! (m1 < m3).all() );
157     VERIFY(! (m1 > m3).all() );
158   }
159   VERIFY(!(m1 > m2 && m1 < m2).any());
160   VERIFY((m1 <= m2 || m1 >= m2).all());
161 
162   // comparisons array to scalar
163   VERIFY( (m1 != (m1(r,c)+1) ).any() );
164   VERIFY( (m1 >  (m1(r,c)-1) ).any() );
165   VERIFY( (m1 <  (m1(r,c)+1) ).any() );
166   VERIFY( (m1 ==  m1(r,c)    ).any() );
167 
168   // comparisons scalar to array
169   VERIFY( ( (m1(r,c)+1) != m1).any() );
170   VERIFY( ( (m1(r,c)-1) <  m1).any() );
171   VERIFY( ( (m1(r,c)+1) >  m1).any() );
172   VERIFY( (  m1(r,c)    == m1).any() );
173 
174   // test Select
175   VERIFY_IS_APPROX( (m1<m2).select(m1,m2), m1.cwiseMin(m2) );
176   VERIFY_IS_APPROX( (m1>m2).select(m1,m2), m1.cwiseMax(m2) );
177   Scalar mid = (m1.cwiseAbs().minCoeff() + m1.cwiseAbs().maxCoeff())/Scalar(2);
178   for (int j=0; j<cols; ++j)
179   for (int i=0; i<rows; ++i)
180     m3(i,j) = abs(m1(i,j))<mid ? 0 : m1(i,j);
181   VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
182                         .select(ArrayType::Zero(rows,cols),m1), m3);
183   // shorter versions:
184   VERIFY_IS_APPROX( (m1.abs()<ArrayType::Constant(rows,cols,mid))
185                         .select(0,m1), m3);
186   VERIFY_IS_APPROX( (m1.abs()>=ArrayType::Constant(rows,cols,mid))
187                         .select(m1,0), m3);
188   // even shorter version:
189   VERIFY_IS_APPROX( (m1.abs()<mid).select(0,m1), m3);
190 
191   // count
192   VERIFY(((m1.abs()+1)>RealScalar(0.1)).count() == rows*cols);
193 
194   // and/or
195   VERIFY( (m1<RealScalar(0) && m1>RealScalar(0)).count() == 0);
196   VERIFY( (m1<RealScalar(0) || m1>=RealScalar(0)).count() == rows*cols);
197   RealScalar a = m1.abs().mean();
198   VERIFY( (m1<-a || m1>a).count() == (m1.abs()>a).count());
199 
200   typedef Array<typename ArrayType::Index, Dynamic, 1> ArrayOfIndices;
201 
202   // TODO allows colwise/rowwise for array
203   VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).colwise().count(), ArrayOfIndices::Constant(cols,rows).transpose());
204   VERIFY_IS_APPROX(((m1.abs()+1)>RealScalar(0.1)).rowwise().count(), ArrayOfIndices::Constant(rows, cols));
205 }
206 
array_real(const ArrayType & m)207 template<typename ArrayType> void array_real(const ArrayType& m)
208 {
209   using std::abs;
210   using std::sqrt;
211   typedef typename ArrayType::Index Index;
212   typedef typename ArrayType::Scalar Scalar;
213   typedef typename NumTraits<Scalar>::Real RealScalar;
214 
215   Index rows = m.rows();
216   Index cols = m.cols();
217 
218   ArrayType m1 = ArrayType::Random(rows, cols),
219             m2 = ArrayType::Random(rows, cols),
220             m3(rows, cols),
221             m4 = m1;
222 
223   m4 = (m4.abs()==Scalar(0)).select(1,m4);
224 
225   Scalar  s1 = internal::random<Scalar>();
226 
227   // these tests are mostly to check possible compilation issues with free-functions.
228   VERIFY_IS_APPROX(m1.sin(), sin(m1));
229   VERIFY_IS_APPROX(m1.cos(), cos(m1));
230   VERIFY_IS_APPROX(m1.tan(), tan(m1));
231   VERIFY_IS_APPROX(m1.asin(), asin(m1));
232   VERIFY_IS_APPROX(m1.acos(), acos(m1));
233   VERIFY_IS_APPROX(m1.atan(), atan(m1));
234   VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
235   VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
236   VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
237 
238   VERIFY_IS_APPROX(m1.arg(), arg(m1));
239   VERIFY_IS_APPROX(m1.round(), round(m1));
240   VERIFY_IS_APPROX(m1.floor(), floor(m1));
241   VERIFY_IS_APPROX(m1.ceil(), ceil(m1));
242   VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
243   VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
244   VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
245   VERIFY_IS_APPROX(m1.inverse(), inverse(m1));
246   VERIFY_IS_APPROX(m1.abs(), abs(m1));
247   VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
248   VERIFY_IS_APPROX(m1.square(), square(m1));
249   VERIFY_IS_APPROX(m1.cube(), cube(m1));
250   VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
251   VERIFY_IS_APPROX(m1.sign(), sign(m1));
252 
253 
254   // avoid NaNs with abs() so verification doesn't fail
255   m3 = m1.abs();
256   VERIFY_IS_APPROX(m3.sqrt(), sqrt(abs(m1)));
257   VERIFY_IS_APPROX(m3.rsqrt(), Scalar(1)/sqrt(abs(m1)));
258   VERIFY_IS_APPROX(rsqrt(m3), Scalar(1)/sqrt(abs(m1)));
259   VERIFY_IS_APPROX(m3.log(), log(m3));
260   VERIFY_IS_APPROX(m3.log1p(), log1p(m3));
261   VERIFY_IS_APPROX(m3.log10(), log10(m3));
262 
263 
264   VERIFY((!(m1>m2) == (m1<=m2)).all());
265 
266   VERIFY_IS_APPROX(sin(m1.asin()), m1);
267   VERIFY_IS_APPROX(cos(m1.acos()), m1);
268   VERIFY_IS_APPROX(tan(m1.atan()), m1);
269   VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
270   VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
271   VERIFY_IS_APPROX(tanh(m1), (0.5*(exp(m1)-exp(-m1)))/(0.5*(exp(m1)+exp(-m1))));
272   VERIFY_IS_APPROX(arg(m1), ((m1<0).template cast<Scalar>())*std::acos(-1.0));
273   VERIFY((round(m1) <= ceil(m1) && round(m1) >= floor(m1)).all());
274   VERIFY((Eigen::isnan)((m1*0.0)/0.0).all());
275   VERIFY((Eigen::isinf)(m4/0.0).all());
276   VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*0.0/0.0)) && (!(Eigen::isfinite)(m4/0.0))).all());
277   VERIFY_IS_APPROX(inverse(inverse(m1)),m1);
278   VERIFY((abs(m1) == m1 || abs(m1) == -m1).all());
279   VERIFY_IS_APPROX(m3, sqrt(abs2(m1)));
280   VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
281   VERIFY_IS_APPROX( m1*m1.sign(),m1.abs());
282   VERIFY_IS_APPROX(m1.sign() * m1.abs(), m1);
283 
284   VERIFY_IS_APPROX(numext::abs2(numext::real(m1)) + numext::abs2(numext::imag(m1)), numext::abs2(m1));
285   VERIFY_IS_APPROX(numext::abs2(real(m1)) + numext::abs2(imag(m1)), numext::abs2(m1));
286   if(!NumTraits<Scalar>::IsComplex)
287     VERIFY_IS_APPROX(numext::real(m1), m1);
288 
289   // shift argument of logarithm so that it is not zero
290   Scalar smallNumber = NumTraits<Scalar>::dummy_precision();
291   VERIFY_IS_APPROX((m3 + smallNumber).log() , log(abs(m1) + smallNumber));
292   VERIFY_IS_APPROX((m3 + smallNumber + 1).log() , log1p(abs(m1) + smallNumber));
293 
294   VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
295   VERIFY_IS_APPROX(m1.exp(), exp(m1));
296   VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
297 
298   VERIFY_IS_APPROX(m3.pow(RealScalar(0.5)), m3.sqrt());
299   VERIFY_IS_APPROX(pow(m3,RealScalar(0.5)), m3.sqrt());
300 
301   VERIFY_IS_APPROX(m3.pow(RealScalar(-0.5)), m3.rsqrt());
302   VERIFY_IS_APPROX(pow(m3,RealScalar(-0.5)), m3.rsqrt());
303 
304   VERIFY_IS_APPROX(log10(m3), log(m3)/log(10));
305 
306   // scalar by array division
307   const RealScalar tiny = sqrt(std::numeric_limits<RealScalar>::epsilon());
308   s1 += Scalar(tiny);
309   m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
310   VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
311 
312   // check inplace transpose
313   m3 = m1;
314   m3.transposeInPlace();
315   VERIFY_IS_APPROX(m3, m1.transpose());
316   m3.transposeInPlace();
317   VERIFY_IS_APPROX(m3, m1);
318 }
319 
array_complex(const ArrayType & m)320 template<typename ArrayType> void array_complex(const ArrayType& m)
321 {
322   typedef typename ArrayType::Index Index;
323   typedef typename ArrayType::Scalar Scalar;
324   typedef typename NumTraits<Scalar>::Real RealScalar;
325 
326   Index rows = m.rows();
327   Index cols = m.cols();
328 
329   ArrayType m1 = ArrayType::Random(rows, cols),
330             m2(rows, cols),
331             m4 = m1;
332 
333   m4.real() = (m4.real().abs()==RealScalar(0)).select(RealScalar(1),m4.real());
334   m4.imag() = (m4.imag().abs()==RealScalar(0)).select(RealScalar(1),m4.imag());
335 
336   Array<RealScalar, -1, -1> m3(rows, cols);
337 
338   for (Index i = 0; i < m.rows(); ++i)
339     for (Index j = 0; j < m.cols(); ++j)
340       m2(i,j) = sqrt(m1(i,j));
341 
342   // these tests are mostly to check possible compilation issues with free-functions.
343   VERIFY_IS_APPROX(m1.sin(), sin(m1));
344   VERIFY_IS_APPROX(m1.cos(), cos(m1));
345   VERIFY_IS_APPROX(m1.tan(), tan(m1));
346   VERIFY_IS_APPROX(m1.sinh(), sinh(m1));
347   VERIFY_IS_APPROX(m1.cosh(), cosh(m1));
348   VERIFY_IS_APPROX(m1.tanh(), tanh(m1));
349   VERIFY_IS_APPROX(m1.arg(), arg(m1));
350   VERIFY((m1.isNaN() == (Eigen::isnan)(m1)).all());
351   VERIFY((m1.isInf() == (Eigen::isinf)(m1)).all());
352   VERIFY((m1.isFinite() == (Eigen::isfinite)(m1)).all());
353   VERIFY_IS_APPROX(m1.inverse(), inverse(m1));
354   VERIFY_IS_APPROX(m1.log(), log(m1));
355   VERIFY_IS_APPROX(m1.log10(), log10(m1));
356   VERIFY_IS_APPROX(m1.abs(), abs(m1));
357   VERIFY_IS_APPROX(m1.abs2(), abs2(m1));
358   VERIFY_IS_APPROX(m1.sqrt(), sqrt(m1));
359   VERIFY_IS_APPROX(m1.square(), square(m1));
360   VERIFY_IS_APPROX(m1.cube(), cube(m1));
361   VERIFY_IS_APPROX(cos(m1+RealScalar(3)*m2), cos((m1+RealScalar(3)*m2).eval()));
362   VERIFY_IS_APPROX(m1.sign(), sign(m1));
363 
364 
365   VERIFY_IS_APPROX(m1.exp() * m2.exp(), exp(m1+m2));
366   VERIFY_IS_APPROX(m1.exp(), exp(m1));
367   VERIFY_IS_APPROX(m1.exp() / m2.exp(),(m1-m2).exp());
368 
369   VERIFY_IS_APPROX(sinh(m1), 0.5*(exp(m1)-exp(-m1)));
370   VERIFY_IS_APPROX(cosh(m1), 0.5*(exp(m1)+exp(-m1)));
371   VERIFY_IS_APPROX(tanh(m1), (0.5*(exp(m1)-exp(-m1)))/(0.5*(exp(m1)+exp(-m1))));
372 
373   for (Index i = 0; i < m.rows(); ++i)
374     for (Index j = 0; j < m.cols(); ++j)
375       m3(i,j) = std::atan2(imag(m1(i,j)), real(m1(i,j)));
376   VERIFY_IS_APPROX(arg(m1), m3);
377 
378   std::complex<RealScalar> zero(0.0,0.0);
379   VERIFY((Eigen::isnan)(m1*zero/zero).all());
380 #if EIGEN_COMP_MSVC
381   // msvc complex division is not robust
382   VERIFY((Eigen::isinf)(m4/RealScalar(0)).all());
383 #else
384 #if EIGEN_COMP_CLANG
385   // clang's complex division is notoriously broken too
386   if((numext::isinf)(m4(0,0)/RealScalar(0))) {
387 #endif
388     VERIFY((Eigen::isinf)(m4/zero).all());
389 #if EIGEN_COMP_CLANG
390   }
391   else
392   {
393     VERIFY((Eigen::isinf)(m4.real()/zero.real()).all());
394   }
395 #endif
396 #endif // MSVC
397 
398   VERIFY(((Eigen::isfinite)(m1) && (!(Eigen::isfinite)(m1*zero/zero)) && (!(Eigen::isfinite)(m1/zero))).all());
399 
400   VERIFY_IS_APPROX(inverse(inverse(m1)),m1);
401   VERIFY_IS_APPROX(conj(m1.conjugate()), m1);
402   VERIFY_IS_APPROX(abs(m1), sqrt(square(real(m1))+square(imag(m1))));
403   VERIFY_IS_APPROX(abs(m1), sqrt(abs2(m1)));
404   VERIFY_IS_APPROX(log10(m1), log(m1)/log(10));
405 
406   VERIFY_IS_APPROX( m1.sign(), -(-m1).sign() );
407   VERIFY_IS_APPROX( m1.sign() * m1.abs(), m1);
408 
409   // scalar by array division
410   Scalar  s1 = internal::random<Scalar>();
411   const RealScalar tiny = std::sqrt(std::numeric_limits<RealScalar>::epsilon());
412   s1 += Scalar(tiny);
413   m1 += ArrayType::Constant(rows,cols,Scalar(tiny));
414   VERIFY_IS_APPROX(s1/m1, s1 * m1.inverse());
415 
416   // check inplace transpose
417   m2 = m1;
418   m2.transposeInPlace();
419   VERIFY_IS_APPROX(m2, m1.transpose());
420   m2.transposeInPlace();
421   VERIFY_IS_APPROX(m2, m1);
422 
423 }
424 
min_max(const ArrayType & m)425 template<typename ArrayType> void min_max(const ArrayType& m)
426 {
427   typedef typename ArrayType::Index Index;
428   typedef typename ArrayType::Scalar Scalar;
429 
430   Index rows = m.rows();
431   Index cols = m.cols();
432 
433   ArrayType m1 = ArrayType::Random(rows, cols);
434 
435   // min/max with array
436   Scalar maxM1 = m1.maxCoeff();
437   Scalar minM1 = m1.minCoeff();
438 
439   VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)(ArrayType::Constant(rows,cols, minM1)));
440   VERIFY_IS_APPROX(m1, (m1.min)(ArrayType::Constant(rows,cols, maxM1)));
441 
442   VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)(ArrayType::Constant(rows,cols, maxM1)));
443   VERIFY_IS_APPROX(m1, (m1.max)(ArrayType::Constant(rows,cols, minM1)));
444 
445   // min/max with scalar input
446   VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, minM1), (m1.min)( minM1));
447   VERIFY_IS_APPROX(m1, (m1.min)( maxM1));
448 
449   VERIFY_IS_APPROX(ArrayType::Constant(rows,cols, maxM1), (m1.max)( maxM1));
450   VERIFY_IS_APPROX(m1, (m1.max)( minM1));
451 
452 }
453 
test_array()454 void test_array()
455 {
456   for(int i = 0; i < g_repeat; i++) {
457     CALL_SUBTEST_1( array(Array<float, 1, 1>()) );
458     CALL_SUBTEST_2( array(Array22f()) );
459     CALL_SUBTEST_3( array(Array44d()) );
460     CALL_SUBTEST_4( array(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
461     CALL_SUBTEST_5( array(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
462     CALL_SUBTEST_6( array(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
463   }
464   for(int i = 0; i < g_repeat; i++) {
465     CALL_SUBTEST_1( comparisons(Array<float, 1, 1>()) );
466     CALL_SUBTEST_2( comparisons(Array22f()) );
467     CALL_SUBTEST_3( comparisons(Array44d()) );
468     CALL_SUBTEST_5( comparisons(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
469     CALL_SUBTEST_6( comparisons(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
470   }
471   for(int i = 0; i < g_repeat; i++) {
472     CALL_SUBTEST_1( min_max(Array<float, 1, 1>()) );
473     CALL_SUBTEST_2( min_max(Array22f()) );
474     CALL_SUBTEST_3( min_max(Array44d()) );
475     CALL_SUBTEST_5( min_max(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
476     CALL_SUBTEST_6( min_max(ArrayXXi(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
477   }
478   for(int i = 0; i < g_repeat; i++) {
479     CALL_SUBTEST_1( array_real(Array<float, 1, 1>()) );
480     CALL_SUBTEST_2( array_real(Array22f()) );
481     CALL_SUBTEST_3( array_real(Array44d()) );
482     CALL_SUBTEST_5( array_real(ArrayXXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
483   }
484   for(int i = 0; i < g_repeat; i++) {
485     CALL_SUBTEST_4( array_complex(ArrayXXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
486   }
487 
488   VERIFY((internal::is_same< internal::global_math_functions_filtering_base<int>::type, int >::value));
489   VERIFY((internal::is_same< internal::global_math_functions_filtering_base<float>::type, float >::value));
490   VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Array2i>::type, ArrayBase<Array2i> >::value));
491   typedef CwiseUnaryOp<internal::scalar_abs_op<double>, ArrayXd > Xpr;
492   VERIFY((internal::is_same< internal::global_math_functions_filtering_base<Xpr>::type,
493                            ArrayBase<Xpr>
494                          >::value));
495 }
496