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