xref: /aosp_15_r20/external/eigen/test/vectorwiseop.cpp (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2011 Benoit Jacob <[email protected]>
5*bf2c3715SXin Li // Copyright (C) 2015 Gael Guennebaud <[email protected]>
6*bf2c3715SXin Li //
7*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
8*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
9*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10*bf2c3715SXin Li 
11*bf2c3715SXin Li #define TEST_ENABLE_TEMPORARY_TRACKING
12*bf2c3715SXin Li #define EIGEN_NO_STATIC_ASSERT
13*bf2c3715SXin Li 
14*bf2c3715SXin Li #include "main.h"
15*bf2c3715SXin Li 
vectorwiseop_array(const ArrayType & m)16*bf2c3715SXin Li template<typename ArrayType> void vectorwiseop_array(const ArrayType& m)
17*bf2c3715SXin Li {
18*bf2c3715SXin Li   typedef typename ArrayType::Scalar Scalar;
19*bf2c3715SXin Li   typedef Array<Scalar, ArrayType::RowsAtCompileTime, 1> ColVectorType;
20*bf2c3715SXin Li   typedef Array<Scalar, 1, ArrayType::ColsAtCompileTime> RowVectorType;
21*bf2c3715SXin Li 
22*bf2c3715SXin Li   Index rows = m.rows();
23*bf2c3715SXin Li   Index cols = m.cols();
24*bf2c3715SXin Li   Index r = internal::random<Index>(0, rows-1),
25*bf2c3715SXin Li         c = internal::random<Index>(0, cols-1);
26*bf2c3715SXin Li 
27*bf2c3715SXin Li   ArrayType m1 = ArrayType::Random(rows, cols),
28*bf2c3715SXin Li             m2(rows, cols),
29*bf2c3715SXin Li             m3(rows, cols);
30*bf2c3715SXin Li 
31*bf2c3715SXin Li   ColVectorType colvec = ColVectorType::Random(rows);
32*bf2c3715SXin Li   RowVectorType rowvec = RowVectorType::Random(cols);
33*bf2c3715SXin Li 
34*bf2c3715SXin Li   // test addition
35*bf2c3715SXin Li 
36*bf2c3715SXin Li   m2 = m1;
37*bf2c3715SXin Li   m2.colwise() += colvec;
38*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.colwise() + colvec);
39*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.col(c), m1.col(c) + colvec);
40*bf2c3715SXin Li 
41*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m2.colwise() += colvec.transpose());
42*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.colwise() + colvec.transpose());
43*bf2c3715SXin Li 
44*bf2c3715SXin Li   m2 = m1;
45*bf2c3715SXin Li   m2.rowwise() += rowvec;
46*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.rowwise() + rowvec);
47*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.row(r), m1.row(r) + rowvec);
48*bf2c3715SXin Li 
49*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m2.rowwise() += rowvec.transpose());
50*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.rowwise() + rowvec.transpose());
51*bf2c3715SXin Li 
52*bf2c3715SXin Li   // test substraction
53*bf2c3715SXin Li 
54*bf2c3715SXin Li   m2 = m1;
55*bf2c3715SXin Li   m2.colwise() -= colvec;
56*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.colwise() - colvec);
57*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.col(c), m1.col(c) - colvec);
58*bf2c3715SXin Li 
59*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m2.colwise() -= colvec.transpose());
60*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.colwise() - colvec.transpose());
61*bf2c3715SXin Li 
62*bf2c3715SXin Li   m2 = m1;
63*bf2c3715SXin Li   m2.rowwise() -= rowvec;
64*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.rowwise() - rowvec);
65*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.row(r), m1.row(r) - rowvec);
66*bf2c3715SXin Li 
67*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m2.rowwise() -= rowvec.transpose());
68*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.rowwise() - rowvec.transpose());
69*bf2c3715SXin Li 
70*bf2c3715SXin Li   // test multiplication
71*bf2c3715SXin Li 
72*bf2c3715SXin Li   m2 = m1;
73*bf2c3715SXin Li   m2.colwise() *= colvec;
74*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.colwise() * colvec);
75*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.col(c), m1.col(c) * colvec);
76*bf2c3715SXin Li 
77*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m2.colwise() *= colvec.transpose());
78*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.colwise() * colvec.transpose());
79*bf2c3715SXin Li 
80*bf2c3715SXin Li   m2 = m1;
81*bf2c3715SXin Li   m2.rowwise() *= rowvec;
82*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.rowwise() * rowvec);
83*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.row(r), m1.row(r) * rowvec);
84*bf2c3715SXin Li 
85*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m2.rowwise() *= rowvec.transpose());
86*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.rowwise() * rowvec.transpose());
87*bf2c3715SXin Li 
88*bf2c3715SXin Li   // test quotient
89*bf2c3715SXin Li 
90*bf2c3715SXin Li   m2 = m1;
91*bf2c3715SXin Li   m2.colwise() /= colvec;
92*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.colwise() / colvec);
93*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.col(c), m1.col(c) / colvec);
94*bf2c3715SXin Li 
95*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m2.colwise() /= colvec.transpose());
96*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.colwise() / colvec.transpose());
97*bf2c3715SXin Li 
98*bf2c3715SXin Li   m2 = m1;
99*bf2c3715SXin Li   m2.rowwise() /= rowvec;
100*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.rowwise() / rowvec);
101*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.row(r), m1.row(r) / rowvec);
102*bf2c3715SXin Li 
103*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m2.rowwise() /= rowvec.transpose());
104*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.rowwise() / rowvec.transpose());
105*bf2c3715SXin Li 
106*bf2c3715SXin Li   m2 = m1;
107*bf2c3715SXin Li   // yes, there might be an aliasing issue there but ".rowwise() /="
108*bf2c3715SXin Li   // is supposed to evaluate " m2.colwise().sum()" into a temporary to avoid
109*bf2c3715SXin Li   // evaluating the reduction multiple times
110*bf2c3715SXin Li   if(ArrayType::RowsAtCompileTime>2 || ArrayType::RowsAtCompileTime==Dynamic)
111*bf2c3715SXin Li   {
112*bf2c3715SXin Li     m2.rowwise() /= m2.colwise().sum();
113*bf2c3715SXin Li     VERIFY_IS_APPROX(m2, m1.rowwise() / m1.colwise().sum());
114*bf2c3715SXin Li   }
115*bf2c3715SXin Li 
116*bf2c3715SXin Li   // all/any
117*bf2c3715SXin Li   Array<bool,Dynamic,Dynamic> mb(rows,cols);
118*bf2c3715SXin Li   mb = (m1.real()<=0.7).colwise().all();
119*bf2c3715SXin Li   VERIFY( (mb.col(c) == (m1.real().col(c)<=0.7).all()).all() );
120*bf2c3715SXin Li   mb = (m1.real()<=0.7).rowwise().all();
121*bf2c3715SXin Li   VERIFY( (mb.row(r) == (m1.real().row(r)<=0.7).all()).all() );
122*bf2c3715SXin Li 
123*bf2c3715SXin Li   mb = (m1.real()>=0.7).colwise().any();
124*bf2c3715SXin Li   VERIFY( (mb.col(c) == (m1.real().col(c)>=0.7).any()).all() );
125*bf2c3715SXin Li   mb = (m1.real()>=0.7).rowwise().any();
126*bf2c3715SXin Li   VERIFY( (mb.row(r) == (m1.real().row(r)>=0.7).any()).all() );
127*bf2c3715SXin Li }
128*bf2c3715SXin Li 
vectorwiseop_matrix(const MatrixType & m)129*bf2c3715SXin Li template<typename MatrixType> void vectorwiseop_matrix(const MatrixType& m)
130*bf2c3715SXin Li {
131*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
132*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
133*bf2c3715SXin Li   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> ColVectorType;
134*bf2c3715SXin Li   typedef Matrix<Scalar, 1, MatrixType::ColsAtCompileTime> RowVectorType;
135*bf2c3715SXin Li   typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, 1> RealColVectorType;
136*bf2c3715SXin Li   typedef Matrix<RealScalar, 1, MatrixType::ColsAtCompileTime> RealRowVectorType;
137*bf2c3715SXin Li   typedef Matrix<Scalar,Dynamic,Dynamic> MatrixX;
138*bf2c3715SXin Li 
139*bf2c3715SXin Li   Index rows = m.rows();
140*bf2c3715SXin Li   Index cols = m.cols();
141*bf2c3715SXin Li   Index r = internal::random<Index>(0, rows-1),
142*bf2c3715SXin Li         c = internal::random<Index>(0, cols-1);
143*bf2c3715SXin Li 
144*bf2c3715SXin Li   MatrixType m1 = MatrixType::Random(rows, cols),
145*bf2c3715SXin Li             m2(rows, cols),
146*bf2c3715SXin Li             m3(rows, cols);
147*bf2c3715SXin Li 
148*bf2c3715SXin Li   ColVectorType colvec = ColVectorType::Random(rows);
149*bf2c3715SXin Li   RowVectorType rowvec = RowVectorType::Random(cols);
150*bf2c3715SXin Li   RealColVectorType rcres;
151*bf2c3715SXin Li   RealRowVectorType rrres;
152*bf2c3715SXin Li 
153*bf2c3715SXin Li   // test broadcast assignment
154*bf2c3715SXin Li   m2 = m1;
155*bf2c3715SXin Li   m2.colwise() = colvec;
156*bf2c3715SXin Li   for(Index j=0; j<cols; ++j)
157*bf2c3715SXin Li     VERIFY_IS_APPROX(m2.col(j), colvec);
158*bf2c3715SXin Li   m2.rowwise() = rowvec;
159*bf2c3715SXin Li   for(Index i=0; i<rows; ++i)
160*bf2c3715SXin Li     VERIFY_IS_APPROX(m2.row(i), rowvec);
161*bf2c3715SXin Li   if(rows>1)
162*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m2.colwise() = colvec.transpose());
163*bf2c3715SXin Li   if(cols>1)
164*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m2.rowwise() = rowvec.transpose());
165*bf2c3715SXin Li 
166*bf2c3715SXin Li   // test addition
167*bf2c3715SXin Li 
168*bf2c3715SXin Li   m2 = m1;
169*bf2c3715SXin Li   m2.colwise() += colvec;
170*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.colwise() + colvec);
171*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.col(c), m1.col(c) + colvec);
172*bf2c3715SXin Li 
173*bf2c3715SXin Li   if(rows>1)
174*bf2c3715SXin Li   {
175*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m2.colwise() += colvec.transpose());
176*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m1.colwise() + colvec.transpose());
177*bf2c3715SXin Li   }
178*bf2c3715SXin Li 
179*bf2c3715SXin Li   m2 = m1;
180*bf2c3715SXin Li   m2.rowwise() += rowvec;
181*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.rowwise() + rowvec);
182*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.row(r), m1.row(r) + rowvec);
183*bf2c3715SXin Li 
184*bf2c3715SXin Li   if(cols>1)
185*bf2c3715SXin Li   {
186*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m2.rowwise() += rowvec.transpose());
187*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m1.rowwise() + rowvec.transpose());
188*bf2c3715SXin Li   }
189*bf2c3715SXin Li 
190*bf2c3715SXin Li   // test substraction
191*bf2c3715SXin Li 
192*bf2c3715SXin Li   m2 = m1;
193*bf2c3715SXin Li   m2.colwise() -= colvec;
194*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.colwise() - colvec);
195*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.col(c), m1.col(c) - colvec);
196*bf2c3715SXin Li 
197*bf2c3715SXin Li   if(rows>1)
198*bf2c3715SXin Li   {
199*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m2.colwise() -= colvec.transpose());
200*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m1.colwise() - colvec.transpose());
201*bf2c3715SXin Li   }
202*bf2c3715SXin Li 
203*bf2c3715SXin Li   m2 = m1;
204*bf2c3715SXin Li   m2.rowwise() -= rowvec;
205*bf2c3715SXin Li   VERIFY_IS_APPROX(m2, m1.rowwise() - rowvec);
206*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.row(r), m1.row(r) - rowvec);
207*bf2c3715SXin Li 
208*bf2c3715SXin Li   if(cols>1)
209*bf2c3715SXin Li   {
210*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m2.rowwise() -= rowvec.transpose());
211*bf2c3715SXin Li     VERIFY_RAISES_ASSERT(m1.rowwise() - rowvec.transpose());
212*bf2c3715SXin Li   }
213*bf2c3715SXin Li 
214*bf2c3715SXin Li   // ------ partial reductions ------
215*bf2c3715SXin Li 
216*bf2c3715SXin Li   #define TEST_PARTIAL_REDUX_BASIC(FUNC,ROW,COL,PREPROCESS) {                          \
217*bf2c3715SXin Li     ROW = m1 PREPROCESS .colwise().FUNC ;                                              \
218*bf2c3715SXin Li     for(Index k=0; k<cols; ++k) VERIFY_IS_APPROX(ROW(k), m1.col(k) PREPROCESS .FUNC ); \
219*bf2c3715SXin Li     COL = m1 PREPROCESS .rowwise().FUNC ;                                              \
220*bf2c3715SXin Li     for(Index k=0; k<rows; ++k) VERIFY_IS_APPROX(COL(k), m1.row(k) PREPROCESS .FUNC ); \
221*bf2c3715SXin Li   }
222*bf2c3715SXin Li 
223*bf2c3715SXin Li   TEST_PARTIAL_REDUX_BASIC(sum(),        rowvec,colvec,EIGEN_EMPTY);
224*bf2c3715SXin Li   TEST_PARTIAL_REDUX_BASIC(prod(),       rowvec,colvec,EIGEN_EMPTY);
225*bf2c3715SXin Li   TEST_PARTIAL_REDUX_BASIC(mean(),       rowvec,colvec,EIGEN_EMPTY);
226*bf2c3715SXin Li   TEST_PARTIAL_REDUX_BASIC(minCoeff(),   rrres, rcres, .real());
227*bf2c3715SXin Li   TEST_PARTIAL_REDUX_BASIC(maxCoeff(),   rrres, rcres, .real());
228*bf2c3715SXin Li   TEST_PARTIAL_REDUX_BASIC(norm(),       rrres, rcres, EIGEN_EMPTY);
229*bf2c3715SXin Li   TEST_PARTIAL_REDUX_BASIC(squaredNorm(),rrres, rcres, EIGEN_EMPTY);
230*bf2c3715SXin Li   TEST_PARTIAL_REDUX_BASIC(redux(internal::scalar_sum_op<Scalar,Scalar>()),rowvec,colvec,EIGEN_EMPTY);
231*bf2c3715SXin Li 
232*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.cwiseAbs().colwise().sum(), m1.colwise().template lpNorm<1>());
233*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.cwiseAbs().rowwise().sum(), m1.rowwise().template lpNorm<1>());
234*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.cwiseAbs().colwise().maxCoeff(), m1.colwise().template lpNorm<Infinity>());
235*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.cwiseAbs().rowwise().maxCoeff(), m1.rowwise().template lpNorm<Infinity>());
236*bf2c3715SXin Li 
237*bf2c3715SXin Li   // regression for bug 1158
238*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.cwiseAbs().colwise().sum().x(), m1.col(0).cwiseAbs().sum());
239*bf2c3715SXin Li 
240*bf2c3715SXin Li   // test normalized
241*bf2c3715SXin Li   m2 = m1.colwise().normalized();
242*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.col(c), m1.col(c).normalized());
243*bf2c3715SXin Li   m2 = m1.rowwise().normalized();
244*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.row(r), m1.row(r).normalized());
245*bf2c3715SXin Li 
246*bf2c3715SXin Li   // test normalize
247*bf2c3715SXin Li   m2 = m1;
248*bf2c3715SXin Li   m2.colwise().normalize();
249*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.col(c), m1.col(c).normalized());
250*bf2c3715SXin Li   m2 = m1;
251*bf2c3715SXin Li   m2.rowwise().normalize();
252*bf2c3715SXin Li   VERIFY_IS_APPROX(m2.row(r), m1.row(r).normalized());
253*bf2c3715SXin Li 
254*bf2c3715SXin Li   // test with partial reduction of products
255*bf2c3715SXin Li   Matrix<Scalar,MatrixType::RowsAtCompileTime,MatrixType::RowsAtCompileTime> m1m1 = m1 * m1.transpose();
256*bf2c3715SXin Li   VERIFY_IS_APPROX( (m1 * m1.transpose()).colwise().sum(), m1m1.colwise().sum());
257*bf2c3715SXin Li   Matrix<Scalar,1,MatrixType::RowsAtCompileTime> tmp(rows);
258*bf2c3715SXin Li   VERIFY_EVALUATION_COUNT( tmp = (m1 * m1.transpose()).colwise().sum(), 1);
259*bf2c3715SXin Li 
260*bf2c3715SXin Li   m2 = m1.rowwise() - (m1.colwise().sum()/RealScalar(m1.rows())).eval();
261*bf2c3715SXin Li   m1 = m1.rowwise() - (m1.colwise().sum()/RealScalar(m1.rows()));
262*bf2c3715SXin Li   VERIFY_IS_APPROX( m1, m2 );
263*bf2c3715SXin Li   VERIFY_EVALUATION_COUNT( m2 = (m1.rowwise() - m1.colwise().sum()/RealScalar(m1.rows())), (MatrixType::RowsAtCompileTime!=1 ? 1 : 0) );
264*bf2c3715SXin Li 
265*bf2c3715SXin Li   // test empty expressions
266*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleCols(0,0).rowwise().sum().eval(), MatrixX::Zero(rows,1));
267*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleRows(0,0).colwise().sum().eval(), MatrixX::Zero(1,cols));
268*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleCols(0,fix<0>).rowwise().sum().eval(), MatrixX::Zero(rows,1));
269*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleRows(0,fix<0>).colwise().sum().eval(), MatrixX::Zero(1,cols));
270*bf2c3715SXin Li 
271*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleCols(0,0).rowwise().prod().eval(), MatrixX::Ones(rows,1));
272*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleRows(0,0).colwise().prod().eval(), MatrixX::Ones(1,cols));
273*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleCols(0,fix<0>).rowwise().prod().eval(), MatrixX::Ones(rows,1));
274*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleRows(0,fix<0>).colwise().prod().eval(), MatrixX::Ones(1,cols));
275*bf2c3715SXin Li 
276*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.matrix().middleCols(0,0).rowwise().squaredNorm().eval(), MatrixX::Zero(rows,1));
277*bf2c3715SXin Li 
278*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.real().middleCols(0,0).rowwise().minCoeff().eval());
279*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(m1.real().middleRows(0,0).colwise().maxCoeff().eval());
280*bf2c3715SXin Li   VERIFY_IS_EQUAL(m1.real().middleRows(0,0).rowwise().maxCoeff().eval().rows(),0);
281*bf2c3715SXin Li   VERIFY_IS_EQUAL(m1.real().middleCols(0,0).colwise().maxCoeff().eval().cols(),0);
282*bf2c3715SXin Li   VERIFY_IS_EQUAL(m1.real().middleRows(0,fix<0>).rowwise().maxCoeff().eval().rows(),0);
283*bf2c3715SXin Li   VERIFY_IS_EQUAL(m1.real().middleCols(0,fix<0>).colwise().maxCoeff().eval().cols(),0);
284*bf2c3715SXin Li }
285*bf2c3715SXin Li 
EIGEN_DECLARE_TEST(vectorwiseop)286*bf2c3715SXin Li EIGEN_DECLARE_TEST(vectorwiseop)
287*bf2c3715SXin Li {
288*bf2c3715SXin Li   CALL_SUBTEST_1( vectorwiseop_array(Array22cd()) );
289*bf2c3715SXin Li   CALL_SUBTEST_2( vectorwiseop_array(Array<double, 3, 2>()) );
290*bf2c3715SXin Li   CALL_SUBTEST_3( vectorwiseop_array(ArrayXXf(3, 4)) );
291*bf2c3715SXin Li   CALL_SUBTEST_4( vectorwiseop_matrix(Matrix4cf()) );
292*bf2c3715SXin Li   CALL_SUBTEST_5( vectorwiseop_matrix(Matrix4f()) );
293*bf2c3715SXin Li   CALL_SUBTEST_5( vectorwiseop_matrix(Vector4f()) );
294*bf2c3715SXin Li   CALL_SUBTEST_5( vectorwiseop_matrix(Matrix<float,4,5>()) );
295*bf2c3715SXin Li   CALL_SUBTEST_6( vectorwiseop_matrix(MatrixXd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
296*bf2c3715SXin Li   CALL_SUBTEST_7( vectorwiseop_matrix(VectorXd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
297*bf2c3715SXin Li   CALL_SUBTEST_7( vectorwiseop_matrix(RowVectorXd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
298*bf2c3715SXin Li }
299