xref: /aosp_15_r20/external/eigen/test/product_extra.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) 2006-2008 Benoit Jacob <[email protected]>
5*bf2c3715SXin Li //
6*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
7*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
8*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9*bf2c3715SXin Li 
10*bf2c3715SXin Li #include "main.h"
11*bf2c3715SXin Li 
product_extra(const MatrixType & m)12*bf2c3715SXin Li template<typename MatrixType> void product_extra(const MatrixType& m)
13*bf2c3715SXin Li {
14*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
15*bf2c3715SXin Li   typedef Matrix<Scalar, 1, Dynamic> RowVectorType;
16*bf2c3715SXin Li   typedef Matrix<Scalar, Dynamic, 1> ColVectorType;
17*bf2c3715SXin Li   typedef Matrix<Scalar, Dynamic, Dynamic,
18*bf2c3715SXin Li                          MatrixType::Flags&RowMajorBit> OtherMajorMatrixType;
19*bf2c3715SXin Li 
20*bf2c3715SXin Li   Index rows = m.rows();
21*bf2c3715SXin Li   Index cols = m.cols();
22*bf2c3715SXin Li 
23*bf2c3715SXin Li   MatrixType m1 = MatrixType::Random(rows, cols),
24*bf2c3715SXin Li              m2 = MatrixType::Random(rows, cols),
25*bf2c3715SXin Li              m3(rows, cols),
26*bf2c3715SXin Li              mzero = MatrixType::Zero(rows, cols),
27*bf2c3715SXin Li              identity = MatrixType::Identity(rows, rows),
28*bf2c3715SXin Li              square = MatrixType::Random(rows, rows),
29*bf2c3715SXin Li              res = MatrixType::Random(rows, rows),
30*bf2c3715SXin Li              square2 = MatrixType::Random(cols, cols),
31*bf2c3715SXin Li              res2 = MatrixType::Random(cols, cols);
32*bf2c3715SXin Li   RowVectorType v1 = RowVectorType::Random(rows), vrres(rows);
33*bf2c3715SXin Li   ColVectorType vc2 = ColVectorType::Random(cols), vcres(cols);
34*bf2c3715SXin Li   OtherMajorMatrixType tm1 = m1;
35*bf2c3715SXin Li 
36*bf2c3715SXin Li   Scalar s1 = internal::random<Scalar>(),
37*bf2c3715SXin Li          s2 = internal::random<Scalar>(),
38*bf2c3715SXin Li          s3 = internal::random<Scalar>();
39*bf2c3715SXin Li 
40*bf2c3715SXin Li   VERIFY_IS_APPROX(m3.noalias() = m1 * m2.adjoint(),                 m1 * m2.adjoint().eval());
41*bf2c3715SXin Li   VERIFY_IS_APPROX(m3.noalias() = m1.adjoint() * square.adjoint(),   m1.adjoint().eval() * square.adjoint().eval());
42*bf2c3715SXin Li   VERIFY_IS_APPROX(m3.noalias() = m1.adjoint() * m2,                 m1.adjoint().eval() * m2);
43*bf2c3715SXin Li   VERIFY_IS_APPROX(m3.noalias() = (s1 * m1.adjoint()) * m2,          (s1 * m1.adjoint()).eval() * m2);
44*bf2c3715SXin Li   VERIFY_IS_APPROX(m3.noalias() = ((s1 * m1).adjoint()) * m2,        (numext::conj(s1) * m1.adjoint()).eval() * m2);
45*bf2c3715SXin Li   VERIFY_IS_APPROX(m3.noalias() = (- m1.adjoint() * s1) * (s3 * m2), (- m1.adjoint()  * s1).eval() * (s3 * m2).eval());
46*bf2c3715SXin Li   VERIFY_IS_APPROX(m3.noalias() = (s2 * m1.adjoint() * s1) * m2,     (s2 * m1.adjoint()  * s1).eval() * m2);
47*bf2c3715SXin Li   VERIFY_IS_APPROX(m3.noalias() = (-m1*s2) * s1*m2.adjoint(),        (-m1*s2).eval() * (s1*m2.adjoint()).eval());
48*bf2c3715SXin Li 
49*bf2c3715SXin Li   // a very tricky case where a scale factor has to be automatically conjugated:
50*bf2c3715SXin Li   VERIFY_IS_APPROX( m1.adjoint() * (s1*m2).conjugate(), (m1.adjoint()).eval() * ((s1*m2).conjugate()).eval());
51*bf2c3715SXin Li 
52*bf2c3715SXin Li 
53*bf2c3715SXin Li   // test all possible conjugate combinations for the four matrix-vector product cases:
54*bf2c3715SXin Li 
55*bf2c3715SXin Li   VERIFY_IS_APPROX((-m1.conjugate() * s2) * (s1 * vc2),
56*bf2c3715SXin Li                    (-m1.conjugate()*s2).eval() * (s1 * vc2).eval());
57*bf2c3715SXin Li   VERIFY_IS_APPROX((-m1 * s2) * (s1 * vc2.conjugate()),
58*bf2c3715SXin Li                    (-m1*s2).eval() * (s1 * vc2.conjugate()).eval());
59*bf2c3715SXin Li   VERIFY_IS_APPROX((-m1.conjugate() * s2) * (s1 * vc2.conjugate()),
60*bf2c3715SXin Li                    (-m1.conjugate()*s2).eval() * (s1 * vc2.conjugate()).eval());
61*bf2c3715SXin Li 
62*bf2c3715SXin Li   VERIFY_IS_APPROX((s1 * vc2.transpose()) * (-m1.adjoint() * s2),
63*bf2c3715SXin Li                    (s1 * vc2.transpose()).eval() * (-m1.adjoint()*s2).eval());
64*bf2c3715SXin Li   VERIFY_IS_APPROX((s1 * vc2.adjoint()) * (-m1.transpose() * s2),
65*bf2c3715SXin Li                    (s1 * vc2.adjoint()).eval() * (-m1.transpose()*s2).eval());
66*bf2c3715SXin Li   VERIFY_IS_APPROX((s1 * vc2.adjoint()) * (-m1.adjoint() * s2),
67*bf2c3715SXin Li                    (s1 * vc2.adjoint()).eval() * (-m1.adjoint()*s2).eval());
68*bf2c3715SXin Li 
69*bf2c3715SXin Li   VERIFY_IS_APPROX((-m1.adjoint() * s2) * (s1 * v1.transpose()),
70*bf2c3715SXin Li                    (-m1.adjoint()*s2).eval() * (s1 * v1.transpose()).eval());
71*bf2c3715SXin Li   VERIFY_IS_APPROX((-m1.transpose() * s2) * (s1 * v1.adjoint()),
72*bf2c3715SXin Li                    (-m1.transpose()*s2).eval() * (s1 * v1.adjoint()).eval());
73*bf2c3715SXin Li   VERIFY_IS_APPROX((-m1.adjoint() * s2) * (s1 * v1.adjoint()),
74*bf2c3715SXin Li                    (-m1.adjoint()*s2).eval() * (s1 * v1.adjoint()).eval());
75*bf2c3715SXin Li 
76*bf2c3715SXin Li   VERIFY_IS_APPROX((s1 * v1) * (-m1.conjugate() * s2),
77*bf2c3715SXin Li                    (s1 * v1).eval() * (-m1.conjugate()*s2).eval());
78*bf2c3715SXin Li   VERIFY_IS_APPROX((s1 * v1.conjugate()) * (-m1 * s2),
79*bf2c3715SXin Li                    (s1 * v1.conjugate()).eval() * (-m1*s2).eval());
80*bf2c3715SXin Li   VERIFY_IS_APPROX((s1 * v1.conjugate()) * (-m1.conjugate() * s2),
81*bf2c3715SXin Li                    (s1 * v1.conjugate()).eval() * (-m1.conjugate()*s2).eval());
82*bf2c3715SXin Li 
83*bf2c3715SXin Li   VERIFY_IS_APPROX((-m1.adjoint() * s2) * (s1 * v1.adjoint()),
84*bf2c3715SXin Li                    (-m1.adjoint()*s2).eval() * (s1 * v1.adjoint()).eval());
85*bf2c3715SXin Li 
86*bf2c3715SXin Li   // test the vector-matrix product with non aligned starts
87*bf2c3715SXin Li   Index i = internal::random<Index>(0,m1.rows()-2);
88*bf2c3715SXin Li   Index j = internal::random<Index>(0,m1.cols()-2);
89*bf2c3715SXin Li   Index r = internal::random<Index>(1,m1.rows()-i);
90*bf2c3715SXin Li   Index c = internal::random<Index>(1,m1.cols()-j);
91*bf2c3715SXin Li   Index i2 = internal::random<Index>(0,m1.rows()-1);
92*bf2c3715SXin Li   Index j2 = internal::random<Index>(0,m1.cols()-1);
93*bf2c3715SXin Li 
94*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.col(j2).adjoint() * m1.block(0,j,m1.rows(),c), m1.col(j2).adjoint().eval() * m1.block(0,j,m1.rows(),c).eval());
95*bf2c3715SXin Li   VERIFY_IS_APPROX(m1.block(i,0,r,m1.cols()) * m1.row(i2).adjoint(), m1.block(i,0,r,m1.cols()).eval() * m1.row(i2).adjoint().eval());
96*bf2c3715SXin Li 
97*bf2c3715SXin Li   // test negative strides
98*bf2c3715SXin Li   {
99*bf2c3715SXin Li     Map<MatrixType,Unaligned,Stride<Dynamic,Dynamic> > map1(&m1(rows-1,cols-1), rows, cols, Stride<Dynamic,Dynamic>(-m1.outerStride(),-1));
100*bf2c3715SXin Li     Map<MatrixType,Unaligned,Stride<Dynamic,Dynamic> > map2(&m2(rows-1,cols-1), rows, cols, Stride<Dynamic,Dynamic>(-m2.outerStride(),-1));
101*bf2c3715SXin Li     Map<RowVectorType,Unaligned,InnerStride<-1> > mapv1(&v1(v1.size()-1), v1.size(), InnerStride<-1>(-1));
102*bf2c3715SXin Li     Map<ColVectorType,Unaligned,InnerStride<-1> > mapvc2(&vc2(vc2.size()-1), vc2.size(), InnerStride<-1>(-1));
103*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(map1), m1.reverse());
104*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(map2), m2.reverse());
105*bf2c3715SXin Li     VERIFY_IS_APPROX(m3.noalias() = MatrixType(map1) * MatrixType(map2).adjoint(), m1.reverse() * m2.reverse().adjoint());
106*bf2c3715SXin Li     VERIFY_IS_APPROX(m3.noalias() = map1 * map2.adjoint(), m1.reverse() * m2.reverse().adjoint());
107*bf2c3715SXin Li     VERIFY_IS_APPROX(map1 * vc2, m1.reverse() * vc2);
108*bf2c3715SXin Li     VERIFY_IS_APPROX(m1 * mapvc2, m1 * mapvc2);
109*bf2c3715SXin Li     VERIFY_IS_APPROX(map1.adjoint() * v1.transpose(), m1.adjoint().reverse() * v1.transpose());
110*bf2c3715SXin Li     VERIFY_IS_APPROX(m1.adjoint() * mapv1.transpose(), m1.adjoint() * v1.reverse().transpose());
111*bf2c3715SXin Li   }
112*bf2c3715SXin Li 
113*bf2c3715SXin Li   // regression test
114*bf2c3715SXin Li   MatrixType tmp = m1 * m1.adjoint() * s1;
115*bf2c3715SXin Li   VERIFY_IS_APPROX(tmp, m1 * m1.adjoint() * s1);
116*bf2c3715SXin Li 
117*bf2c3715SXin Li   // regression test for bug 1343, assignment to arrays
118*bf2c3715SXin Li   Array<Scalar,Dynamic,1> a1 = m1 * vc2;
119*bf2c3715SXin Li   VERIFY_IS_APPROX(a1.matrix(),m1*vc2);
120*bf2c3715SXin Li   Array<Scalar,Dynamic,1> a2 = s1 * (m1 * vc2);
121*bf2c3715SXin Li   VERIFY_IS_APPROX(a2.matrix(),s1*m1*vc2);
122*bf2c3715SXin Li   Array<Scalar,1,Dynamic> a3 = v1 * m1;
123*bf2c3715SXin Li   VERIFY_IS_APPROX(a3.matrix(),v1*m1);
124*bf2c3715SXin Li   Array<Scalar,Dynamic,Dynamic> a4 = m1 * m2.adjoint();
125*bf2c3715SXin Li   VERIFY_IS_APPROX(a4.matrix(),m1*m2.adjoint());
126*bf2c3715SXin Li }
127*bf2c3715SXin Li 
128*bf2c3715SXin Li // Regression test for bug reported at http://forum.kde.org/viewtopic.php?f=74&t=96947
mat_mat_scalar_scalar_product()129*bf2c3715SXin Li void mat_mat_scalar_scalar_product()
130*bf2c3715SXin Li {
131*bf2c3715SXin Li   Eigen::Matrix2Xd dNdxy(2, 3);
132*bf2c3715SXin Li   dNdxy << -0.5, 0.5, 0,
133*bf2c3715SXin Li            -0.3, 0, 0.3;
134*bf2c3715SXin Li   double det = 6.0, wt = 0.5;
135*bf2c3715SXin Li   VERIFY_IS_APPROX(dNdxy.transpose()*dNdxy*det*wt, det*wt*dNdxy.transpose()*dNdxy);
136*bf2c3715SXin Li }
137*bf2c3715SXin Li 
138*bf2c3715SXin Li template <typename MatrixType>
zero_sized_objects(const MatrixType & m)139*bf2c3715SXin Li void zero_sized_objects(const MatrixType& m)
140*bf2c3715SXin Li {
141*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
142*bf2c3715SXin Li   const int PacketSize  = internal::packet_traits<Scalar>::size;
143*bf2c3715SXin Li   const int PacketSize1 = PacketSize>1 ?  PacketSize-1 : 1;
144*bf2c3715SXin Li   Index rows = m.rows();
145*bf2c3715SXin Li   Index cols = m.cols();
146*bf2c3715SXin Li 
147*bf2c3715SXin Li   {
148*bf2c3715SXin Li     MatrixType res, a(rows,0), b(0,cols);
149*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a*b), MatrixType::Zero(rows,cols) );
150*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a*a.transpose()), MatrixType::Zero(rows,rows) );
151*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=b.transpose()*b), MatrixType::Zero(cols,cols) );
152*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=b.transpose()*a.transpose()), MatrixType::Zero(cols,rows) );
153*bf2c3715SXin Li   }
154*bf2c3715SXin Li 
155*bf2c3715SXin Li   {
156*bf2c3715SXin Li     MatrixType res, a(rows,cols), b(cols,0);
157*bf2c3715SXin Li     res = a*b;
158*bf2c3715SXin Li     VERIFY(res.rows()==rows && res.cols()==0);
159*bf2c3715SXin Li     b.resize(0,rows);
160*bf2c3715SXin Li     res = b*a;
161*bf2c3715SXin Li     VERIFY(res.rows()==0 && res.cols()==cols);
162*bf2c3715SXin Li   }
163*bf2c3715SXin Li 
164*bf2c3715SXin Li   {
165*bf2c3715SXin Li     Matrix<Scalar,PacketSize,0> a;
166*bf2c3715SXin Li     Matrix<Scalar,0,1> b;
167*bf2c3715SXin Li     Matrix<Scalar,PacketSize,1> res;
168*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a*b), MatrixType::Zero(PacketSize,1) );
169*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a.lazyProduct(b)), MatrixType::Zero(PacketSize,1) );
170*bf2c3715SXin Li   }
171*bf2c3715SXin Li 
172*bf2c3715SXin Li   {
173*bf2c3715SXin Li     Matrix<Scalar,PacketSize1,0> a;
174*bf2c3715SXin Li     Matrix<Scalar,0,1> b;
175*bf2c3715SXin Li     Matrix<Scalar,PacketSize1,1> res;
176*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a*b), MatrixType::Zero(PacketSize1,1) );
177*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a.lazyProduct(b)), MatrixType::Zero(PacketSize1,1) );
178*bf2c3715SXin Li   }
179*bf2c3715SXin Li 
180*bf2c3715SXin Li   {
181*bf2c3715SXin Li     Matrix<Scalar,PacketSize,Dynamic> a(PacketSize,0);
182*bf2c3715SXin Li     Matrix<Scalar,Dynamic,1> b(0,1);
183*bf2c3715SXin Li     Matrix<Scalar,PacketSize,1> res;
184*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a*b), MatrixType::Zero(PacketSize,1) );
185*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a.lazyProduct(b)), MatrixType::Zero(PacketSize,1) );
186*bf2c3715SXin Li   }
187*bf2c3715SXin Li 
188*bf2c3715SXin Li   {
189*bf2c3715SXin Li     Matrix<Scalar,PacketSize1,Dynamic> a(PacketSize1,0);
190*bf2c3715SXin Li     Matrix<Scalar,Dynamic,1> b(0,1);
191*bf2c3715SXin Li     Matrix<Scalar,PacketSize1,1> res;
192*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a*b), MatrixType::Zero(PacketSize1,1) );
193*bf2c3715SXin Li     VERIFY_IS_APPROX( (res=a.lazyProduct(b)), MatrixType::Zero(PacketSize1,1) );
194*bf2c3715SXin Li   }
195*bf2c3715SXin Li }
196*bf2c3715SXin Li 
197*bf2c3715SXin Li template<int>
bug_127()198*bf2c3715SXin Li void bug_127()
199*bf2c3715SXin Li {
200*bf2c3715SXin Li   // Bug 127
201*bf2c3715SXin Li   //
202*bf2c3715SXin Li   // a product of the form lhs*rhs with
203*bf2c3715SXin Li   //
204*bf2c3715SXin Li   // lhs:
205*bf2c3715SXin Li   // rows = 1, cols = 4
206*bf2c3715SXin Li   // RowsAtCompileTime = 1, ColsAtCompileTime = -1
207*bf2c3715SXin Li   // MaxRowsAtCompileTime = 1, MaxColsAtCompileTime = 5
208*bf2c3715SXin Li   //
209*bf2c3715SXin Li   // rhs:
210*bf2c3715SXin Li   // rows = 4, cols = 0
211*bf2c3715SXin Li   // RowsAtCompileTime = -1, ColsAtCompileTime = -1
212*bf2c3715SXin Li   // MaxRowsAtCompileTime = 5, MaxColsAtCompileTime = 1
213*bf2c3715SXin Li   //
214*bf2c3715SXin Li   // was failing on a runtime assertion, because it had been mis-compiled as a dot product because Product.h was using the
215*bf2c3715SXin Li   // max-sizes to detect size 1 indicating vectors, and that didn't account for 0-sized object with max-size 1.
216*bf2c3715SXin Li 
217*bf2c3715SXin Li   Matrix<float,1,Dynamic,RowMajor,1,5> a(1,4);
218*bf2c3715SXin Li   Matrix<float,Dynamic,Dynamic,ColMajor,5,1> b(4,0);
219*bf2c3715SXin Li   a*b;
220*bf2c3715SXin Li }
221*bf2c3715SXin Li 
bug_817()222*bf2c3715SXin Li template<int> void bug_817()
223*bf2c3715SXin Li {
224*bf2c3715SXin Li   ArrayXXf B = ArrayXXf::Random(10,10), C;
225*bf2c3715SXin Li   VectorXf x = VectorXf::Random(10);
226*bf2c3715SXin Li   C = (x.transpose()*B.matrix());
227*bf2c3715SXin Li   B = (x.transpose()*B.matrix());
228*bf2c3715SXin Li   VERIFY_IS_APPROX(B,C);
229*bf2c3715SXin Li }
230*bf2c3715SXin Li 
231*bf2c3715SXin Li template<int>
unaligned_objects()232*bf2c3715SXin Li void unaligned_objects()
233*bf2c3715SXin Li {
234*bf2c3715SXin Li   // Regression test for the bug reported here:
235*bf2c3715SXin Li   // http://forum.kde.org/viewtopic.php?f=74&t=107541
236*bf2c3715SXin Li   // Recall the matrix*vector kernel avoid unaligned loads by loading two packets and then reassemble then.
237*bf2c3715SXin Li   // There was a mistake in the computation of the valid range for fully unaligned objects: in some rare cases,
238*bf2c3715SXin Li   // memory was read outside the allocated matrix memory. Though the values were not used, this might raise segfault.
239*bf2c3715SXin Li   for(int m=450;m<460;++m)
240*bf2c3715SXin Li   {
241*bf2c3715SXin Li     for(int n=8;n<12;++n)
242*bf2c3715SXin Li     {
243*bf2c3715SXin Li       MatrixXf M(m, n);
244*bf2c3715SXin Li       VectorXf v1(n), r1(500);
245*bf2c3715SXin Li       RowVectorXf v2(m), r2(16);
246*bf2c3715SXin Li 
247*bf2c3715SXin Li       M.setRandom();
248*bf2c3715SXin Li       v1.setRandom();
249*bf2c3715SXin Li       v2.setRandom();
250*bf2c3715SXin Li       for(int o=0; o<4; ++o)
251*bf2c3715SXin Li       {
252*bf2c3715SXin Li         r1.segment(o,m).noalias() = M * v1;
253*bf2c3715SXin Li         VERIFY_IS_APPROX(r1.segment(o,m), M * MatrixXf(v1));
254*bf2c3715SXin Li         r2.segment(o,n).noalias() = v2 * M;
255*bf2c3715SXin Li         VERIFY_IS_APPROX(r2.segment(o,n), MatrixXf(v2) * M);
256*bf2c3715SXin Li       }
257*bf2c3715SXin Li     }
258*bf2c3715SXin Li   }
259*bf2c3715SXin Li }
260*bf2c3715SXin Li 
261*bf2c3715SXin Li template<typename T>
262*bf2c3715SXin Li EIGEN_DONT_INLINE
test_compute_block_size(Index m,Index n,Index k)263*bf2c3715SXin Li Index test_compute_block_size(Index m, Index n, Index k)
264*bf2c3715SXin Li {
265*bf2c3715SXin Li   Index mc(m), nc(n), kc(k);
266*bf2c3715SXin Li   internal::computeProductBlockingSizes<T,T>(kc, mc, nc);
267*bf2c3715SXin Li   return kc+mc+nc;
268*bf2c3715SXin Li }
269*bf2c3715SXin Li 
270*bf2c3715SXin Li template<typename T>
compute_block_size()271*bf2c3715SXin Li Index compute_block_size()
272*bf2c3715SXin Li {
273*bf2c3715SXin Li   Index ret = 0;
274*bf2c3715SXin Li   ret += test_compute_block_size<T>(0,1,1);
275*bf2c3715SXin Li   ret += test_compute_block_size<T>(1,0,1);
276*bf2c3715SXin Li   ret += test_compute_block_size<T>(1,1,0);
277*bf2c3715SXin Li   ret += test_compute_block_size<T>(0,0,1);
278*bf2c3715SXin Li   ret += test_compute_block_size<T>(0,1,0);
279*bf2c3715SXin Li   ret += test_compute_block_size<T>(1,0,0);
280*bf2c3715SXin Li   ret += test_compute_block_size<T>(0,0,0);
281*bf2c3715SXin Li   return ret;
282*bf2c3715SXin Li }
283*bf2c3715SXin Li 
284*bf2c3715SXin Li template<typename>
aliasing_with_resize()285*bf2c3715SXin Li void aliasing_with_resize()
286*bf2c3715SXin Li {
287*bf2c3715SXin Li   Index m = internal::random<Index>(10,50);
288*bf2c3715SXin Li   Index n = internal::random<Index>(10,50);
289*bf2c3715SXin Li   MatrixXd A, B, C(m,n), D(m,m);
290*bf2c3715SXin Li   VectorXd a, b, c(n);
291*bf2c3715SXin Li   C.setRandom();
292*bf2c3715SXin Li   D.setRandom();
293*bf2c3715SXin Li   c.setRandom();
294*bf2c3715SXin Li   double s = internal::random<double>(1,10);
295*bf2c3715SXin Li 
296*bf2c3715SXin Li   A = C;
297*bf2c3715SXin Li   B = A * A.transpose();
298*bf2c3715SXin Li   A = A * A.transpose();
299*bf2c3715SXin Li   VERIFY_IS_APPROX(A,B);
300*bf2c3715SXin Li 
301*bf2c3715SXin Li   A = C;
302*bf2c3715SXin Li   B = (A * A.transpose())/s;
303*bf2c3715SXin Li   A = (A * A.transpose())/s;
304*bf2c3715SXin Li   VERIFY_IS_APPROX(A,B);
305*bf2c3715SXin Li 
306*bf2c3715SXin Li   A = C;
307*bf2c3715SXin Li   B = (A * A.transpose()) + D;
308*bf2c3715SXin Li   A = (A * A.transpose()) + D;
309*bf2c3715SXin Li   VERIFY_IS_APPROX(A,B);
310*bf2c3715SXin Li 
311*bf2c3715SXin Li   A = C;
312*bf2c3715SXin Li   B = D + (A * A.transpose());
313*bf2c3715SXin Li   A = D + (A * A.transpose());
314*bf2c3715SXin Li   VERIFY_IS_APPROX(A,B);
315*bf2c3715SXin Li 
316*bf2c3715SXin Li   A = C;
317*bf2c3715SXin Li   B = s * (A * A.transpose());
318*bf2c3715SXin Li   A = s * (A * A.transpose());
319*bf2c3715SXin Li   VERIFY_IS_APPROX(A,B);
320*bf2c3715SXin Li 
321*bf2c3715SXin Li   A = C;
322*bf2c3715SXin Li   a = c;
323*bf2c3715SXin Li   b = (A * a)/s;
324*bf2c3715SXin Li   a = (A * a)/s;
325*bf2c3715SXin Li   VERIFY_IS_APPROX(a,b);
326*bf2c3715SXin Li }
327*bf2c3715SXin Li 
328*bf2c3715SXin Li template<int>
bug_1308()329*bf2c3715SXin Li void bug_1308()
330*bf2c3715SXin Li {
331*bf2c3715SXin Li   int n = 10;
332*bf2c3715SXin Li   MatrixXd r(n,n);
333*bf2c3715SXin Li   VectorXd v = VectorXd::Random(n);
334*bf2c3715SXin Li   r = v * RowVectorXd::Ones(n);
335*bf2c3715SXin Li   VERIFY_IS_APPROX(r, v.rowwise().replicate(n));
336*bf2c3715SXin Li   r = VectorXd::Ones(n) * v.transpose();
337*bf2c3715SXin Li   VERIFY_IS_APPROX(r, v.rowwise().replicate(n).transpose());
338*bf2c3715SXin Li 
339*bf2c3715SXin Li   Matrix4d ones44 = Matrix4d::Ones();
340*bf2c3715SXin Li   Matrix4d m44 = Matrix4d::Ones() * Matrix4d::Ones();
341*bf2c3715SXin Li   VERIFY_IS_APPROX(m44,Matrix4d::Constant(4));
342*bf2c3715SXin Li   VERIFY_IS_APPROX(m44.noalias()=ones44*Matrix4d::Ones(), Matrix4d::Constant(4));
343*bf2c3715SXin Li   VERIFY_IS_APPROX(m44.noalias()=ones44.transpose()*Matrix4d::Ones(), Matrix4d::Constant(4));
344*bf2c3715SXin Li   VERIFY_IS_APPROX(m44.noalias()=Matrix4d::Ones()*ones44, Matrix4d::Constant(4));
345*bf2c3715SXin Li   VERIFY_IS_APPROX(m44.noalias()=Matrix4d::Ones()*ones44.transpose(), Matrix4d::Constant(4));
346*bf2c3715SXin Li 
347*bf2c3715SXin Li   typedef Matrix<double,4,4,RowMajor> RMatrix4d;
348*bf2c3715SXin Li   RMatrix4d r44 = Matrix4d::Ones() * Matrix4d::Ones();
349*bf2c3715SXin Li   VERIFY_IS_APPROX(r44,Matrix4d::Constant(4));
350*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias()=ones44*Matrix4d::Ones(), Matrix4d::Constant(4));
351*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias()=ones44.transpose()*Matrix4d::Ones(), Matrix4d::Constant(4));
352*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias()=Matrix4d::Ones()*ones44, Matrix4d::Constant(4));
353*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias()=Matrix4d::Ones()*ones44.transpose(), Matrix4d::Constant(4));
354*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias()=ones44*RMatrix4d::Ones(), Matrix4d::Constant(4));
355*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias()=ones44.transpose()*RMatrix4d::Ones(), Matrix4d::Constant(4));
356*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias()=RMatrix4d::Ones()*ones44, Matrix4d::Constant(4));
357*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias()=RMatrix4d::Ones()*ones44.transpose(), Matrix4d::Constant(4));
358*bf2c3715SXin Li 
359*bf2c3715SXin Li //   RowVector4d r4;
360*bf2c3715SXin Li   m44.setOnes();
361*bf2c3715SXin Li   r44.setZero();
362*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias() += m44.row(0).transpose() * RowVector4d::Ones(), ones44);
363*bf2c3715SXin Li   r44.setZero();
364*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias() += m44.col(0) * RowVector4d::Ones(), ones44);
365*bf2c3715SXin Li   r44.setZero();
366*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias() += Vector4d::Ones() * m44.row(0), ones44);
367*bf2c3715SXin Li   r44.setZero();
368*bf2c3715SXin Li   VERIFY_IS_APPROX(r44.noalias() += Vector4d::Ones() * m44.col(0).transpose(), ones44);
369*bf2c3715SXin Li }
370*bf2c3715SXin Li 
EIGEN_DECLARE_TEST(product_extra)371*bf2c3715SXin Li EIGEN_DECLARE_TEST(product_extra)
372*bf2c3715SXin Li {
373*bf2c3715SXin Li   for(int i = 0; i < g_repeat; i++) {
374*bf2c3715SXin Li     CALL_SUBTEST_1( product_extra(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
375*bf2c3715SXin Li     CALL_SUBTEST_2( product_extra(MatrixXd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
376*bf2c3715SXin Li     CALL_SUBTEST_2( mat_mat_scalar_scalar_product() );
377*bf2c3715SXin Li     CALL_SUBTEST_3( product_extra(MatrixXcf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2), internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) );
378*bf2c3715SXin Li     CALL_SUBTEST_4( product_extra(MatrixXcd(internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2), internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2))) );
379*bf2c3715SXin Li     CALL_SUBTEST_1( zero_sized_objects(MatrixXf(internal::random<int>(1,EIGEN_TEST_MAX_SIZE), internal::random<int>(1,EIGEN_TEST_MAX_SIZE))) );
380*bf2c3715SXin Li   }
381*bf2c3715SXin Li   CALL_SUBTEST_5( bug_127<0>() );
382*bf2c3715SXin Li   CALL_SUBTEST_5( bug_817<0>() );
383*bf2c3715SXin Li   CALL_SUBTEST_5( bug_1308<0>() );
384*bf2c3715SXin Li   CALL_SUBTEST_6( unaligned_objects<0>() );
385*bf2c3715SXin Li   CALL_SUBTEST_7( compute_block_size<float>() );
386*bf2c3715SXin Li   CALL_SUBTEST_7( compute_block_size<double>() );
387*bf2c3715SXin Li   CALL_SUBTEST_7( compute_block_size<std::complex<double> >() );
388*bf2c3715SXin Li   CALL_SUBTEST_8( aliasing_with_resize<void>() );
389*bf2c3715SXin Li 
390*bf2c3715SXin Li }
391