1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2006-2008 Benoit Jacob <jacob.benoit.1@gmail.com>
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 #include <Eigen/QR>
12
13 template<typename Derived1, typename Derived2>
14 bool areNotApprox(const MatrixBase<Derived1>& m1, const MatrixBase<Derived2>& m2, typename Derived1::RealScalar epsilon = NumTraits<typename Derived1::RealScalar>::dummy_precision())
15 {
16 return !((m1-m2).cwiseAbs2().maxCoeff() < epsilon * epsilon
17 * (std::max)(m1.cwiseAbs2().maxCoeff(), m2.cwiseAbs2().maxCoeff()));
18 }
19
product(const MatrixType & m)20 template<typename MatrixType> void product(const MatrixType& m)
21 {
22 /* this test covers the following files:
23 Identity.h Product.h
24 */
25 typedef typename MatrixType::Scalar Scalar;
26 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> RowVectorType;
27 typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, 1> ColVectorType;
28 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> RowSquareMatrixType;
29 typedef Matrix<Scalar, MatrixType::ColsAtCompileTime, MatrixType::ColsAtCompileTime> ColSquareMatrixType;
30 typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::ColsAtCompileTime,
31 MatrixType::Flags&RowMajorBit?ColMajor:RowMajor> OtherMajorMatrixType;
32
33 Index rows = m.rows();
34 Index cols = m.cols();
35
36 // this test relies a lot on Random.h, and there's not much more that we can do
37 // to test it, hence I consider that we will have tested Random.h
38 MatrixType m1 = MatrixType::Random(rows, cols),
39 m2 = MatrixType::Random(rows, cols),
40 m3(rows, cols);
41 RowSquareMatrixType
42 identity = RowSquareMatrixType::Identity(rows, rows),
43 square = RowSquareMatrixType::Random(rows, rows),
44 res = RowSquareMatrixType::Random(rows, rows);
45 ColSquareMatrixType
46 square2 = ColSquareMatrixType::Random(cols, cols),
47 res2 = ColSquareMatrixType::Random(cols, cols);
48 RowVectorType v1 = RowVectorType::Random(rows);
49 ColVectorType vc2 = ColVectorType::Random(cols), vcres(cols);
50 OtherMajorMatrixType tm1 = m1;
51
52 Scalar s1 = internal::random<Scalar>();
53
54 Index r = internal::random<Index>(0, rows-1),
55 c = internal::random<Index>(0, cols-1),
56 c2 = internal::random<Index>(0, cols-1);
57
58 // begin testing Product.h: only associativity for now
59 // (we use Transpose.h but this doesn't count as a test for it)
60 VERIFY_IS_APPROX((m1*m1.transpose())*m2, m1*(m1.transpose()*m2));
61 m3 = m1;
62 m3 *= m1.transpose() * m2;
63 VERIFY_IS_APPROX(m3, m1 * (m1.transpose()*m2));
64 VERIFY_IS_APPROX(m3, m1 * (m1.transpose()*m2));
65
66 // continue testing Product.h: distributivity
67 VERIFY_IS_APPROX(square*(m1 + m2), square*m1+square*m2);
68 VERIFY_IS_APPROX(square*(m1 - m2), square*m1-square*m2);
69
70 // continue testing Product.h: compatibility with ScalarMultiple.h
71 VERIFY_IS_APPROX(s1*(square*m1), (s1*square)*m1);
72 VERIFY_IS_APPROX(s1*(square*m1), square*(m1*s1));
73
74 // test Product.h together with Identity.h
75 VERIFY_IS_APPROX(v1, identity*v1);
76 VERIFY_IS_APPROX(v1.transpose(), v1.transpose() * identity);
77 // again, test operator() to check const-qualification
78 VERIFY_IS_APPROX(MatrixType::Identity(rows, cols)(r,c), static_cast<Scalar>(r==c));
79
80 if (rows!=cols)
81 VERIFY_RAISES_ASSERT(m3 = m1*m1);
82
83 // test the previous tests were not screwed up because operator* returns 0
84 // (we use the more accurate default epsilon)
85 if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
86 {
87 VERIFY(areNotApprox(m1.transpose()*m2,m2.transpose()*m1));
88 }
89
90 // test optimized operator+= path
91 res = square;
92 res.noalias() += m1 * m2.transpose();
93 VERIFY_IS_APPROX(res, square + m1 * m2.transpose());
94 if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
95 {
96 VERIFY(areNotApprox(res,square + m2 * m1.transpose()));
97 }
98 vcres = vc2;
99 vcres.noalias() += m1.transpose() * v1;
100 VERIFY_IS_APPROX(vcres, vc2 + m1.transpose() * v1);
101
102 // test optimized operator-= path
103 res = square;
104 res.noalias() -= m1 * m2.transpose();
105 VERIFY_IS_APPROX(res, square - (m1 * m2.transpose()));
106 if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
107 {
108 VERIFY(areNotApprox(res,square - m2 * m1.transpose()));
109 }
110 vcres = vc2;
111 vcres.noalias() -= m1.transpose() * v1;
112 VERIFY_IS_APPROX(vcres, vc2 - m1.transpose() * v1);
113
114 // test d ?= a+b*c rules
115 res.noalias() = square + m1 * m2.transpose();
116 VERIFY_IS_APPROX(res, square + m1 * m2.transpose());
117 res.noalias() += square + m1 * m2.transpose();
118 VERIFY_IS_APPROX(res, 2*(square + m1 * m2.transpose()));
119 res.noalias() -= square + m1 * m2.transpose();
120 VERIFY_IS_APPROX(res, square + m1 * m2.transpose());
121
122 // test d ?= a-b*c rules
123 res.noalias() = square - m1 * m2.transpose();
124 VERIFY_IS_APPROX(res, square - m1 * m2.transpose());
125 res.noalias() += square - m1 * m2.transpose();
126 VERIFY_IS_APPROX(res, 2*(square - m1 * m2.transpose()));
127 res.noalias() -= square - m1 * m2.transpose();
128 VERIFY_IS_APPROX(res, square - m1 * m2.transpose());
129
130
131 tm1 = m1;
132 VERIFY_IS_APPROX(tm1.transpose() * v1, m1.transpose() * v1);
133 VERIFY_IS_APPROX(v1.transpose() * tm1, v1.transpose() * m1);
134
135 // test submatrix and matrix/vector product
136 for (int i=0; i<rows; ++i)
137 res.row(i) = m1.row(i) * m2.transpose();
138 VERIFY_IS_APPROX(res, m1 * m2.transpose());
139 // the other way round:
140 for (int i=0; i<rows; ++i)
141 res.col(i) = m1 * m2.transpose().col(i);
142 VERIFY_IS_APPROX(res, m1 * m2.transpose());
143
144 res2 = square2;
145 res2.noalias() += m1.transpose() * m2;
146 VERIFY_IS_APPROX(res2, square2 + m1.transpose() * m2);
147 if (!NumTraits<Scalar>::IsInteger && (std::min)(rows,cols)>1)
148 {
149 VERIFY(areNotApprox(res2,square2 + m2.transpose() * m1));
150 }
151
152 VERIFY_IS_APPROX(res.col(r).noalias() = square.adjoint() * square.col(r), (square.adjoint() * square.col(r)).eval());
153 VERIFY_IS_APPROX(res.col(r).noalias() = square * square.col(r), (square * square.col(r)).eval());
154
155 // vector at runtime (see bug 1166)
156 {
157 RowSquareMatrixType ref(square);
158 ColSquareMatrixType ref2(square2);
159 ref = res = square;
160 VERIFY_IS_APPROX(res.block(0,0,1,rows).noalias() = m1.col(0).transpose() * square.transpose(), (ref.row(0) = m1.col(0).transpose() * square.transpose()));
161 VERIFY_IS_APPROX(res.block(0,0,1,rows).noalias() = m1.block(0,0,rows,1).transpose() * square.transpose(), (ref.row(0) = m1.col(0).transpose() * square.transpose()));
162 VERIFY_IS_APPROX(res.block(0,0,1,rows).noalias() = m1.col(0).transpose() * square, (ref.row(0) = m1.col(0).transpose() * square));
163 VERIFY_IS_APPROX(res.block(0,0,1,rows).noalias() = m1.block(0,0,rows,1).transpose() * square, (ref.row(0) = m1.col(0).transpose() * square));
164 ref2 = res2 = square2;
165 VERIFY_IS_APPROX(res2.block(0,0,1,cols).noalias() = m1.row(0) * square2.transpose(), (ref2.row(0) = m1.row(0) * square2.transpose()));
166 VERIFY_IS_APPROX(res2.block(0,0,1,cols).noalias() = m1.block(0,0,1,cols) * square2.transpose(), (ref2.row(0) = m1.row(0) * square2.transpose()));
167 VERIFY_IS_APPROX(res2.block(0,0,1,cols).noalias() = m1.row(0) * square2, (ref2.row(0) = m1.row(0) * square2));
168 VERIFY_IS_APPROX(res2.block(0,0,1,cols).noalias() = m1.block(0,0,1,cols) * square2, (ref2.row(0) = m1.row(0) * square2));
169 }
170
171 // vector.block() (see bug 1283)
172 {
173 RowVectorType w1(rows);
174 VERIFY_IS_APPROX(square * v1.block(0,0,rows,1), square * v1);
175 VERIFY_IS_APPROX(w1.noalias() = square * v1.block(0,0,rows,1), square * v1);
176 VERIFY_IS_APPROX(w1.block(0,0,rows,1).noalias() = square * v1.block(0,0,rows,1), square * v1);
177
178 Matrix<Scalar,1,MatrixType::ColsAtCompileTime> w2(cols);
179 VERIFY_IS_APPROX(vc2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
180 VERIFY_IS_APPROX(w2.noalias() = vc2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
181 VERIFY_IS_APPROX(w2.block(0,0,1,cols).noalias() = vc2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
182
183 vc2 = square2.block(0,0,1,cols).transpose();
184 VERIFY_IS_APPROX(square2.block(0,0,1,cols) * square2, vc2.transpose() * square2);
185 VERIFY_IS_APPROX(w2.noalias() = square2.block(0,0,1,cols) * square2, vc2.transpose() * square2);
186 VERIFY_IS_APPROX(w2.block(0,0,1,cols).noalias() = square2.block(0,0,1,cols) * square2, vc2.transpose() * square2);
187
188 vc2 = square2.block(0,0,cols,1);
189 VERIFY_IS_APPROX(square2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
190 VERIFY_IS_APPROX(w2.noalias() = square2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
191 VERIFY_IS_APPROX(w2.block(0,0,1,cols).noalias() = square2.block(0,0,cols,1).transpose() * square2, vc2.transpose() * square2);
192 }
193
194 // inner product
195 {
196 Scalar x = square2.row(c) * square2.col(c2);
197 VERIFY_IS_APPROX(x, square2.row(c).transpose().cwiseProduct(square2.col(c2)).sum());
198 }
199
200 // outer product
201 {
202 VERIFY_IS_APPROX(m1.col(c) * m1.row(r), m1.block(0,c,rows,1) * m1.block(r,0,1,cols));
203 VERIFY_IS_APPROX(m1.row(r).transpose() * m1.col(c).transpose(), m1.block(r,0,1,cols).transpose() * m1.block(0,c,rows,1).transpose());
204 VERIFY_IS_APPROX(m1.block(0,c,rows,1) * m1.row(r), m1.block(0,c,rows,1) * m1.block(r,0,1,cols));
205 VERIFY_IS_APPROX(m1.col(c) * m1.block(r,0,1,cols), m1.block(0,c,rows,1) * m1.block(r,0,1,cols));
206 VERIFY_IS_APPROX(m1.leftCols(1) * m1.row(r), m1.block(0,0,rows,1) * m1.block(r,0,1,cols));
207 VERIFY_IS_APPROX(m1.col(c) * m1.topRows(1), m1.block(0,c,rows,1) * m1.block(0,0,1,cols));
208 }
209
210 // Aliasing
211 {
212 ColVectorType x(cols); x.setRandom();
213 ColVectorType z(x);
214 ColVectorType y(cols); y.setZero();
215 ColSquareMatrixType A(cols,cols); A.setRandom();
216 // CwiseBinaryOp
217 VERIFY_IS_APPROX(x = y + A*x, A*z);
218 x = z;
219 // CwiseUnaryOp
220 VERIFY_IS_APPROX(x = Scalar(1.)*(A*x), A*z);
221 }
222
223 // regression for blas_trais
224 {
225 VERIFY_IS_APPROX(square * (square*square).transpose(), square * square.transpose() * square.transpose());
226 VERIFY_IS_APPROX(square * (-(square*square)), -square * square * square);
227 VERIFY_IS_APPROX(square * (s1*(square*square)), s1 * square * square * square);
228 VERIFY_IS_APPROX(square * (square*square).conjugate(), square * square.conjugate() * square.conjugate());
229 }
230
231 }
232