blob: 9d559f5bfea0dd3e3018906808d166c2e3bb9777 [file] [log] [blame]
Gael Guennebaud709e9032009-01-07 17:01:57 +00001// This file is part of Eigen, a lightweight C++ template library
Benoit Jacob6347b1d2009-05-22 20:25:33 +02002// for linear algebra.
Gael Guennebaud709e9032009-01-07 17:01:57 +00003//
Gael Guennebaud22c76092011-03-22 11:58:22 +01004// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
Gael Guennebaud709e9032009-01-07 17:01:57 +00005//
Benoit Jacob69124cf2012-07-13 14:42:47 -04006// 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/.
Gael Guennebaud709e9032009-01-07 17:01:57 +00009
10#include "sparse.h"
11
12template<typename Scalar> void sparse_vector(int rows, int cols)
13{
Gael Guennebaud42e25782011-08-19 14:18:05 +020014 double densityMat = (std::max)(8./(rows*cols), 0.01);
15 double densityVec = (std::max)(8./float(rows), 0.1);
Gael Guennebaud709e9032009-01-07 17:01:57 +000016 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
17 typedef Matrix<Scalar,Dynamic,1> DenseVector;
18 typedef SparseVector<Scalar> SparseVectorType;
19 typedef SparseMatrix<Scalar> SparseMatrixType;
20 Scalar eps = 1e-6;
21
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010022 SparseMatrixType m1(rows,rows);
Gael Guennebaud709e9032009-01-07 17:01:57 +000023 SparseVectorType v1(rows), v2(rows), v3(rows);
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010024 DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
Gael Guennebaud709e9032009-01-07 17:01:57 +000025 DenseVector refV1 = DenseVector::Random(rows),
26 refV2 = DenseVector::Random(rows),
27 refV3 = DenseVector::Random(rows);
28
29 std::vector<int> zerocoords, nonzerocoords;
30 initSparse<Scalar>(densityVec, refV1, v1, &zerocoords, &nonzerocoords);
31 initSparse<Scalar>(densityMat, refM1, m1);
32
33 initSparse<Scalar>(densityVec, refV2, v2);
34 initSparse<Scalar>(densityVec, refV3, v3);
35
Benoit Jacob47160402010-10-25 10:15:22 -040036 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud709e9032009-01-07 17:01:57 +000037
38 // test coeff and coeffRef
39 for (unsigned int i=0; i<zerocoords.size(); ++i)
40 {
41 VERIFY_IS_MUCH_SMALLER_THAN( v1.coeff(zerocoords[i]), eps );
Gael Guennebaud36c478c2009-01-19 22:29:28 +000042 //VERIFY_RAISES_ASSERT( v1.coeffRef(zerocoords[i]) = 5 );
Gael Guennebaud709e9032009-01-07 17:01:57 +000043 }
44 {
45 VERIFY(int(nonzerocoords.size()) == v1.nonZeros());
46 int j=0;
47 for (typename SparseVectorType::InnerIterator it(v1); it; ++it,++j)
48 {
49 VERIFY(nonzerocoords[j]==it.index());
Gael Guennebaudc4c70662009-01-14 14:24:10 +000050 VERIFY(it.value()==v1.coeff(it.index()));
51 VERIFY(it.value()==refV1.coeff(it.index()));
Gael Guennebaud709e9032009-01-07 17:01:57 +000052 }
53 }
54 VERIFY_IS_APPROX(v1, refV1);
55
56 v1.coeffRef(nonzerocoords[0]) = Scalar(5);
57 refV1.coeffRef(nonzerocoords[0]) = Scalar(5);
58 VERIFY_IS_APPROX(v1, refV1);
59
60 VERIFY_IS_APPROX(v1+v2, refV1+refV2);
61 VERIFY_IS_APPROX(v1+v2+v3, refV1+refV2+refV3);
62
63 VERIFY_IS_APPROX(v1*s1-v2, refV1*s1-refV2);
64
Gael Guennebaud2d534662009-01-14 21:27:54 +000065 VERIFY_IS_APPROX(v1*=s1, refV1*=s1);
66 VERIFY_IS_APPROX(v1/=s1, refV1/=s1);
Gael Guennebauda2324d62010-04-13 10:40:55 +020067
Gael Guennebaude7c48fa2009-01-23 13:59:32 +000068 VERIFY_IS_APPROX(v1+=v2, refV1+=refV2);
69 VERIFY_IS_APPROX(v1-=v2, refV1-=refV2);
Gael Guennebaud2d534662009-01-14 21:27:54 +000070
Gael Guennebaud709e9032009-01-07 17:01:57 +000071 VERIFY_IS_APPROX(v1.dot(v2), refV1.dot(refV2));
Gael Guennebauda9688f02009-02-09 09:59:30 +000072 VERIFY_IS_APPROX(v1.dot(refV2), refV1.dot(refV2));
Gael Guennebaud709e9032009-01-07 17:01:57 +000073
Gael Guennebaud2c03e6f2011-12-22 14:01:06 +010074 VERIFY_IS_APPROX(v1.dot(m1*v2), refV1.dot(refM1*refV2));
75 int i = internal::random<int>(0,rows-1);
76 VERIFY_IS_APPROX(v1.dot(m1.col(i)), refV1.dot(refM1.col(i)));
77
78
Gael Guennebauda2324d62010-04-13 10:40:55 +020079 VERIFY_IS_APPROX(v1.squaredNorm(), refV1.squaredNorm());
80
Gael Guennebaude75b1eb2012-07-25 09:33:50 +020081 // test aliasing
82 VERIFY_IS_APPROX((v1 = -v1), (refV1 = -refV1));
83 VERIFY_IS_APPROX((v1 = v1.transpose()), (refV1 = refV1.transpose().eval()));
84 VERIFY_IS_APPROX((v1 += -v1), (refV1 += -refV1));
85
Gael Guennebaud709e9032009-01-07 17:01:57 +000086}
87
88void test_sparse_vector()
89{
90 for(int i = 0; i < g_repeat; i++) {
Benoit Jacob2840ac72009-10-28 18:19:29 -040091 CALL_SUBTEST_1( sparse_vector<double>(8, 8) );
92 CALL_SUBTEST_2( sparse_vector<std::complex<double> >(16, 16) );
93 CALL_SUBTEST_1( sparse_vector<double>(299, 535) );
Gael Guennebaud709e9032009-01-07 17:01:57 +000094 }
95}
96