blob: 6a3e19105a3e5e175e94f56ff586ce7ebdcb1048 [file] [log] [blame]
Gael Guennebaud86ccd992008-11-05 13:47:55 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaud86ccd992008-11-05 13:47:55 +00003//
Gael Guennebaud22c76092011-03-22 11:58:22 +01004// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
Gael Guennebaud86ccd992008-11-05 13:47:55 +00005// Copyright (C) 2008 Daniel Gomez Ferro <dgomezferro@gmail.com>
6//
7// Eigen is free software; you can redistribute it and/or
8// modify it under the terms of the GNU Lesser General Public
9// License as published by the Free Software Foundation; either
10// version 3 of the License, or (at your option) any later version.
11//
12// Alternatively, you can redistribute it and/or
13// modify it under the terms of the GNU General Public License as
14// published by the Free Software Foundation; either version 2 of
15// the License, or (at your option) any later version.
16//
17// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
18// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
19// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
20// GNU General Public License for more details.
21//
22// You should have received a copy of the GNU Lesser General Public
23// License and a copy of the GNU General Public License along with
24// Eigen. If not, see <http://www.gnu.org/licenses/>.
25
26#include "sparse.h"
27
Gael Guennebaud178858f2009-01-19 15:20:45 +000028template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref)
29{
Hauke Heibelf1679c72010-06-20 17:37:56 +020030 typedef typename SparseMatrixType::Index Index;
31
32 const Index rows = ref.rows();
33 const Index cols = ref.cols();
Gael Guennebaud178858f2009-01-19 15:20:45 +000034 typedef typename SparseMatrixType::Scalar Scalar;
35 enum { Flags = SparseMatrixType::Flags };
Gael Guennebaud9f795582009-12-16 19:18:40 +010036
Gael Guennebaud42e25782011-08-19 14:18:05 +020037 double density = (std::max)(8./(rows*cols), 0.01);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000038 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
39 typedef Matrix<Scalar,Dynamic,1> DenseVector;
40 Scalar eps = 1e-6;
41
Gael Guennebaud178858f2009-01-19 15:20:45 +000042 SparseMatrixType m(rows, cols);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000043 DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
44 DenseVector vec1 = DenseVector::Random(rows);
Benoit Jacob47160402010-10-25 10:15:22 -040045 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud86ccd992008-11-05 13:47:55 +000046
47 std::vector<Vector2i> zeroCoords;
48 std::vector<Vector2i> nonzeroCoords;
49 initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
Gael Guennebaud9f795582009-12-16 19:18:40 +010050
Gael Guennebaud86ccd992008-11-05 13:47:55 +000051 if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
52 return;
53
54 // test coeff and coeffRef
55 for (int i=0; i<(int)zeroCoords.size(); ++i)
56 {
57 VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
Hauke Heibel7bc8e3a2010-10-25 22:13:49 +020058 if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
Gael Guennebaud178858f2009-01-19 15:20:45 +000059 VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
Gael Guennebaud86ccd992008-11-05 13:47:55 +000060 }
61 VERIFY_IS_APPROX(m, refMat);
62
63 m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
64 refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
65
66 VERIFY_IS_APPROX(m, refMat);
Gael Guennebaudc4c70662009-01-14 14:24:10 +000067 /*
Gael Guennebaud86ccd992008-11-05 13:47:55 +000068 // test InnerIterators and Block expressions
69 for (int t=0; t<10; ++t)
70 {
Benoit Jacob47160402010-10-25 10:15:22 -040071 int j = internal::random<int>(0,cols-1);
72 int i = internal::random<int>(0,rows-1);
73 int w = internal::random<int>(1,cols-j-1);
74 int h = internal::random<int>(1,rows-i-1);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000075
Gael Guennebaudc4c70662009-01-14 14:24:10 +000076// VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
Gael Guennebaud86ccd992008-11-05 13:47:55 +000077 for(int c=0; c<w; c++)
78 {
79 VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
80 for(int r=0; r<h; r++)
81 {
Gael Guennebaudc4c70662009-01-14 14:24:10 +000082// VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r));
Gael Guennebaud86ccd992008-11-05 13:47:55 +000083 }
84 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +000085// for(int r=0; r<h; r++)
86// {
87// VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
88// for(int c=0; c<w; c++)
89// {
90// VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
91// }
92// }
Gael Guennebaud86ccd992008-11-05 13:47:55 +000093 }
94
95 for(int c=0; c<cols; c++)
96 {
97 VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
98 VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
99 }
100
101 for(int r=0; r<rows; r++)
102 {
103 VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
104 VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
105 }
106 */
107
Gael Guennebaud28293142009-05-04 14:25:12 +0000108 // test insert (inner random)
Gael Guennebaud5015e482008-12-11 18:26:24 +0000109 {
110 DenseMatrix m1(rows,cols);
111 m1.setZero();
Gael Guennebaud178858f2009-01-19 15:20:45 +0000112 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100113 if(internal::random<int>()%2)
114 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud5015e482008-12-11 18:26:24 +0000115 for (int j=0; j<cols; ++j)
116 {
117 for (int k=0; k<rows/2; ++k)
118 {
Benoit Jacob47160402010-10-25 10:15:22 -0400119 int i = internal::random<int>(0,rows-1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000120 if (m1.coeff(i,j)==Scalar(0))
Benoit Jacob47160402010-10-25 10:15:22 -0400121 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud5015e482008-12-11 18:26:24 +0000122 }
123 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000124 m2.finalize();
125 VERIFY_IS_APPROX(m2,m1);
126 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100127
Gael Guennebaud28293142009-05-04 14:25:12 +0000128 // test insert (fully random)
129 {
130 DenseMatrix m1(rows,cols);
131 m1.setZero();
132 SparseMatrixType m2(rows,cols);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100133 if(internal::random<int>()%2)
134 m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
Gael Guennebaud28293142009-05-04 14:25:12 +0000135 for (int k=0; k<rows*cols; ++k)
136 {
Benoit Jacob47160402010-10-25 10:15:22 -0400137 int i = internal::random<int>(0,rows-1);
138 int j = internal::random<int>(0,cols-1);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100139 if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
Benoit Jacob47160402010-10-25 10:15:22 -0400140 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100141 else
142 {
143 Scalar v = internal::random<Scalar>();
144 m2.coeffRef(i,j) += v;
145 m1(i,j) += v;
146 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000147 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000148 VERIFY_IS_APPROX(m2,m1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000149 }
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200150
151 // test insert (un-compressed)
152 for(int mode=0;mode<4;++mode)
153 {
154 DenseMatrix m1(rows,cols);
155 m1.setZero();
156 SparseMatrixType m2(rows,cols);
157 VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
158 m2.reserve(r);
159 for (int k=0; k<rows*cols; ++k)
160 {
161 int i = internal::random<int>(0,rows-1);
162 int j = internal::random<int>(0,cols-1);
163 if (m1.coeff(i,j)==Scalar(0))
164 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
165 if(mode==3)
166 m2.reserve(r);
167 }
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100168 if(internal::random<int>()%2)
169 m2.makeCompressed();
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200170 VERIFY_IS_APPROX(m2,m1);
171 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100172
Gael Guennebaud2d534662009-01-14 21:27:54 +0000173 // test basic computations
174 {
175 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
176 DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
177 DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
178 DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud178858f2009-01-19 15:20:45 +0000179 SparseMatrixType m1(rows, rows);
180 SparseMatrixType m2(rows, rows);
181 SparseMatrixType m3(rows, rows);
182 SparseMatrixType m4(rows, rows);
Gael Guennebaud2d534662009-01-14 21:27:54 +0000183 initSparse<Scalar>(density, refM1, m1);
184 initSparse<Scalar>(density, refM2, m2);
185 initSparse<Scalar>(density, refM3, m3);
186 initSparse<Scalar>(density, refM4, m4);
187
188 VERIFY_IS_APPROX(m1+m2, refM1+refM2);
189 VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
Gael Guennebaud9f795582009-12-16 19:18:40 +0100190 VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
Gael Guennebaud2d534662009-01-14 21:27:54 +0000191 VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
192
193 VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
194 VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
Gael Guennebaud9f795582009-12-16 19:18:40 +0100195
Gael Guennebaude7c48fa2009-01-23 13:59:32 +0000196 VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
197 VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
Gael Guennebaud9f795582009-12-16 19:18:40 +0100198
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100199 if(SparseMatrixType::IsRowMajor)
200 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
201 else
202 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
Gael Guennebaud9f795582009-12-16 19:18:40 +0100203
Gael Guennebaud91e392a2011-12-03 23:49:37 +0100204 VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
205 VERIFY_IS_APPROX(m1.real(), refM1.real());
206
Gael Guennebaud2d534662009-01-14 21:27:54 +0000207 refM4.setRandom();
208 // sparse cwise* dense
Gael Guennebaud9f795582009-12-16 19:18:40 +0100209 VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
Gael Guennebaud2d534662009-01-14 21:27:54 +0000210// VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
211 }
212
Gael Guennebaudece48a62010-06-18 11:28:30 +0200213 // test transpose
214 {
215 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
216 SparseMatrixType m2(rows, rows);
217 initSparse<Scalar>(density, refMat2, m2);
218 VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
219 VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
220
221 VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
222 }
223
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000224 // test innerVector()
225 {
226 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud178858f2009-01-19 15:20:45 +0000227 SparseMatrixType m2(rows, rows);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000228 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200229 int j0 = internal::random<int>(0,rows-1);
230 int j1 = internal::random<int>(0,rows-1);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100231 if(SparseMatrixType::IsRowMajor)
232 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
233 else
234 VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
235
236 if(SparseMatrixType::IsRowMajor)
237 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
238 else
239 VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100240
241 SparseMatrixType m3(rows,rows);
242 m3.reserve(VectorXi::Constant(rows,rows/2));
243 for(int j=0; j<rows; ++j)
244 for(int k=0; k<j; ++k)
245 m3.insertByOuterInner(j,k) = k+1;
246 for(int j=0; j<rows; ++j)
247 {
Gael Guennebaud67ae94f2011-12-23 09:41:14 +0100248 VERIFY(j==internal::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100249 if(j>0)
Gael Guennebaud67ae94f2011-12-23 09:41:14 +0100250 VERIFY(j==internal::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100251 }
252 m3.makeCompressed();
253 for(int j=0; j<rows; ++j)
254 {
Gael Guennebaud67ae94f2011-12-23 09:41:14 +0100255 VERIFY(j==internal::real(m3.innerVector(j).nonZeros()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100256 if(j>0)
Gael Guennebaud67ae94f2011-12-23 09:41:14 +0100257 VERIFY(j==internal::real(m3.innerVector(j).lastCoeff()));
Gael Guennebaud50d756b2011-12-20 18:10:02 +0100258 }
259
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000260 //m2.innerVector(j0) = 2*m2.innerVector(j1);
261 //refMat2.col(j0) = 2*refMat2.col(j1);
262 //VERIFY_IS_APPROX(m2, refMat2);
263 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100264
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000265 // test innerVectors()
266 {
267 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
268 SparseMatrixType m2(rows, rows);
269 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200270 int j0 = internal::random<int>(0,rows-2);
271 int j1 = internal::random<int>(0,rows-2);
Gael Guennebaud42e25782011-08-19 14:18:05 +0200272 int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100273 if(SparseMatrixType::IsRowMajor)
274 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
275 else
276 VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
277 if(SparseMatrixType::IsRowMajor)
278 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
279 refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
280 else
281 VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
282 refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
Gael Guennebaud8ce45032009-01-27 22:48:17 +0000283 //m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
284 //refMat2.block(0,j0,rows,n0) = refMat2.block(0,j0,rows,n0) + refMat2.block(0,j1,rows,n0);
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000285 }
286
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000287 // test prune
288 {
289 SparseMatrixType m2(rows, rows);
290 DenseMatrix refM2(rows, rows);
291 refM2.setZero();
292 int countFalseNonZero = 0;
293 int countTrueNonZero = 0;
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000294 for (int j=0; j<m2.outerSize(); ++j)
Gael Guennebaud28293142009-05-04 14:25:12 +0000295 {
296 m2.startVec(j);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000297 for (int i=0; i<m2.innerSize(); ++i)
298 {
Benoit Jacob47160402010-10-25 10:15:22 -0400299 float x = internal::random<float>(0,1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000300 if (x<0.1)
301 {
302 // do nothing
303 }
304 else if (x<0.5)
305 {
306 countFalseNonZero++;
Gael Guennebaud8710bd22010-06-02 13:32:13 +0200307 m2.insertBackByOuterInner(j,i) = Scalar(0);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000308 }
309 else
310 {
311 countTrueNonZero++;
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100312 m2.insertBackByOuterInner(j,i) = Scalar(1);
313 if(SparseMatrixType::IsRowMajor)
314 refM2(j,i) = Scalar(1);
315 else
316 refM2(i,j) = Scalar(1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000317 }
318 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000319 }
320 m2.finalize();
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000321 VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
322 VERIFY_IS_APPROX(m2, refM2);
Hauke Heibeld204ec42010-11-02 14:33:33 +0100323 m2.prune(Scalar(1));
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000324 VERIFY(countTrueNonZero==m2.nonZeros());
325 VERIFY_IS_APPROX(m2, refM2);
326 }
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100327
328 // test triangularView
329 {
330 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
331 SparseMatrixType m2(rows, rows), m3(rows, rows);
332 initSparse<Scalar>(density, refMat2, m2);
333 refMat3 = refMat2.template triangularView<Lower>();
334 m3 = m2.template triangularView<Lower>();
335 VERIFY_IS_APPROX(m3, refMat3);
336
337 refMat3 = refMat2.template triangularView<Upper>();
338 m3 = m2.template triangularView<Upper>();
339 VERIFY_IS_APPROX(m3, refMat3);
340
341 refMat3 = refMat2.template triangularView<UnitUpper>();
342 m3 = m2.template triangularView<UnitUpper>();
343 VERIFY_IS_APPROX(m3, refMat3);
344
345 refMat3 = refMat2.template triangularView<UnitLower>();
346 m3 = m2.template triangularView<UnitLower>();
347 VERIFY_IS_APPROX(m3, refMat3);
348 }
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200349
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100350 // test selfadjointView
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100351 if(!SparseMatrixType::IsRowMajor)
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100352 {
353 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
354 SparseMatrixType m2(rows, rows), m3(rows, rows);
355 initSparse<Scalar>(density, refMat2, m2);
356 refMat3 = refMat2.template selfadjointView<Lower>();
357 m3 = m2.template selfadjointView<Lower>();
358 VERIFY_IS_APPROX(m3, refMat3);
359 }
360
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200361 // test sparseView
362 {
363 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
364 SparseMatrixType m2(rows, rows);
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100365 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200366 VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
367 }
Gael Guennebaud82f9aa12011-12-04 21:49:21 +0100368
369 // test diagonal
370 {
371 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
372 SparseMatrixType m2(rows, rows);
373 initSparse<Scalar>(density, refMat2, m2);
374 VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
375 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000376}
377
378void test_sparse_basic()
379{
380 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200381 int s = Eigen::internal::random<int>(1,50);
382 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100383 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(s, s)) ));
384 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(s, s)) ));
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200385 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
386 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100387 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,long int>(s, s)) ));
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000388 }
389}