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