blob: 17897fd618cb52296047b79376f080b6982482bc [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>
Desire NUENTSA4cd82452013-06-11 14:42:29 +02006// Copyright (C) 2013 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
Gael Guennebaud86ccd992008-11-05 13:47:55 +00007//
Benoit Jacob69124cf2012-07-13 14:42:47 -04008// This Source Code Form is subject to the terms of the Mozilla
9// Public License v. 2.0. If a copy of the MPL was not distributed
10// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
Gael Guennebaud86ccd992008-11-05 13:47:55 +000011
Gael Guennebaud8214cf12018-10-11 10:27:23 +020012#ifndef EIGEN_SPARSE_TEST_INCLUDED_FROM_SPARSE_EXTRA
Gael Guennebaudc43154b2015-03-04 10:16:46 +010013static long g_realloc_count = 0;
14#define EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN g_realloc_count++;
15
Gael Guennebaudeec0dfd2018-10-10 22:50:15 +020016static long g_dense_op_sparse_count = 0;
17#define EIGEN_SPARSE_ASSIGNMENT_FROM_DENSE_OP_SPARSE_PLUGIN g_dense_op_sparse_count++;
18#define EIGEN_SPARSE_ASSIGNMENT_FROM_SPARSE_ADD_DENSE_PLUGIN g_dense_op_sparse_count+=10;
19#define EIGEN_SPARSE_ASSIGNMENT_FROM_SPARSE_SUB_DENSE_PLUGIN g_dense_op_sparse_count+=20;
Gael Guennebaud8214cf12018-10-11 10:27:23 +020020#endif
Gael Guennebaudeec0dfd2018-10-10 22:50:15 +020021
Gael Guennebaud86ccd992008-11-05 13:47:55 +000022#include "sparse.h"
23
Gael Guennebaud178858f2009-01-19 15:20:45 +000024template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref)
25{
Christoph Hertzberge8cdbed2014-12-04 22:48:53 +010026 typedef typename SparseMatrixType::StorageIndex StorageIndex;
27 typedef Matrix<StorageIndex,2,1> Vector2;
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +020028
Gael Guennebaudfc202ba2015-02-13 18:57:41 +010029 const Index rows = ref.rows();
30 const Index cols = ref.cols();
Gael Guennebaud0a537cb2016-02-12 15:58:31 +010031 //const Index inner = ref.innerSize();
32 //const Index outer = ref.outerSize();
Christoph Hertzberg0833b822014-10-31 17:12:13 +010033
Gael Guennebaud178858f2009-01-19 15:20:45 +000034 typedef typename SparseMatrixType::Scalar Scalar;
Gael Guennebaud71362672016-12-27 16:34:30 +010035 typedef typename SparseMatrixType::RealScalar RealScalar;
Gael Guennebaud178858f2009-01-19 15:20:45 +000036 enum { Flags = SparseMatrixType::Flags };
Gael Guennebaud9f795582009-12-16 19:18:40 +010037
Gael Guennebaud42e25782011-08-19 14:18:05 +020038 double density = (std::max)(8./(rows*cols), 0.01);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000039 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
40 typedef Matrix<Scalar,Dynamic,1> DenseVector;
41 Scalar eps = 1e-6;
42
Benoit Jacob47160402010-10-25 10:15:22 -040043 Scalar s1 = internal::random<Scalar>();
Gael Guennebaud86ccd992008-11-05 13:47:55 +000044 {
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020045 SparseMatrixType m(rows, cols);
46 DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
47 DenseVector vec1 = DenseVector::Random(rows);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000048
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +020049 std::vector<Vector2> zeroCoords;
50 std::vector<Vector2> nonzeroCoords;
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020051 initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000052
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020053 // test coeff and coeffRef
Christoph Hertzberg0833b822014-10-31 17:12:13 +010054 for (std::size_t i=0; i<zeroCoords.size(); ++i)
Gael Guennebaud86ccd992008-11-05 13:47:55 +000055 {
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020056 VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
57 if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
Christoph Hertzberg0833b822014-10-31 17:12:13 +010058 VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[i].x(),zeroCoords[i].y()) = 5 );
Gael Guennebaud86ccd992008-11-05 13:47:55 +000059 }
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020060 VERIFY_IS_APPROX(m, refMat);
Gael Guennebaud86ccd992008-11-05 13:47:55 +000061
Christoph Hertzberg0833b822014-10-31 17:12:13 +010062 if(!nonzeroCoords.empty()) {
63 m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
64 refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
65 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +000066
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020067 VERIFY_IS_APPROX(m, refMat);
Christoph Hertzberg0833b822014-10-31 17:12:13 +010068
Gael Guennebauda915f022013-06-28 16:16:02 +020069 // test assertion
70 VERIFY_RAISES_ASSERT( m.coeffRef(-1,1) = 0 );
71 VERIFY_RAISES_ASSERT( m.coeffRef(0,m.cols()) = 0 );
Gael Guennebaudd1d7a1a2013-06-23 19:11:32 +020072 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +000073
Gael Guennebaud28293142009-05-04 14:25:12 +000074 // test insert (inner random)
Gael Guennebaud5015e482008-12-11 18:26:24 +000075 {
76 DenseMatrix m1(rows,cols);
77 m1.setZero();
Gael Guennebaud178858f2009-01-19 15:20:45 +000078 SparseMatrixType m2(rows,cols);
Gael Guennebaudc43154b2015-03-04 10:16:46 +010079 bool call_reserve = internal::random<int>()%2;
80 Index nnz = internal::random<int>(1,int(rows)/2);
81 if(call_reserve)
82 {
83 if(internal::random<int>()%2)
84 m2.reserve(VectorXi::Constant(m2.outerSize(), int(nnz)));
85 else
86 m2.reserve(m2.outerSize() * nnz);
87 }
88 g_realloc_count = 0;
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +020089 for (Index j=0; j<cols; ++j)
Gael Guennebaud5015e482008-12-11 18:26:24 +000090 {
Gael Guennebaudc43154b2015-03-04 10:16:46 +010091 for (Index k=0; k<nnz; ++k)
Gael Guennebaud5015e482008-12-11 18:26:24 +000092 {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +020093 Index i = internal::random<Index>(0,rows-1);
Gael Guennebaud5015e482008-12-11 18:26:24 +000094 if (m1.coeff(i,j)==Scalar(0))
Benoit Jacob47160402010-10-25 10:15:22 -040095 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud5015e482008-12-11 18:26:24 +000096 }
97 }
Gael Guennebaudc43154b2015-03-04 10:16:46 +010098
99 if(call_reserve && !SparseMatrixType::IsRowMajor)
100 {
101 VERIFY(g_realloc_count==0);
102 }
103
Gael Guennebaud28293142009-05-04 14:25:12 +0000104 m2.finalize();
105 VERIFY_IS_APPROX(m2,m1);
106 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100107
Gael Guennebaud28293142009-05-04 14:25:12 +0000108 // test insert (fully random)
109 {
110 DenseMatrix m1(rows,cols);
111 m1.setZero();
112 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 Guennebaud28293142009-05-04 14:25:12 +0000115 for (int k=0; k<rows*cols; ++k)
116 {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200117 Index i = internal::random<Index>(0,rows-1);
118 Index j = internal::random<Index>(0,cols-1);
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100119 if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
Benoit Jacob47160402010-10-25 10:15:22 -0400120 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100121 else
122 {
123 Scalar v = internal::random<Scalar>();
124 m2.coeffRef(i,j) += v;
125 m1(i,j) += v;
126 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000127 }
Gael Guennebaudc4c70662009-01-14 14:24:10 +0000128 VERIFY_IS_APPROX(m2,m1);
Gael Guennebaud5015e482008-12-11 18:26:24 +0000129 }
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200130
131 // test insert (un-compressed)
132 for(int mode=0;mode<4;++mode)
133 {
134 DenseMatrix m1(rows,cols);
135 m1.setZero();
136 SparseMatrixType m2(rows,cols);
Gael Guennebaudb47ef142014-07-08 16:47:11 +0200137 VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? int(m2.innerSize()) : std::max<int>(1,int(m2.innerSize())/8)));
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200138 m2.reserve(r);
Christoph Hertzberge8cdbed2014-12-04 22:48:53 +0100139 for (Index k=0; k<rows*cols; ++k)
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200140 {
Gael Guennebaud6d1f5db2013-07-10 23:48:26 +0200141 Index i = internal::random<Index>(0,rows-1);
142 Index j = internal::random<Index>(0,cols-1);
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200143 if (m1.coeff(i,j)==Scalar(0))
144 m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
145 if(mode==3)
146 m2.reserve(r);
147 }
Gael Guennebaud4ca89f32011-12-02 19:00:16 +0100148 if(internal::random<int>()%2)
149 m2.makeCompressed();
Gael Guennebaud7706baf2011-09-08 13:42:54 +0200150 VERIFY_IS_APPROX(m2,m1);
151 }
Gael Guennebaud9f795582009-12-16 19:18:40 +0100152
Gael Guennebaud4e602832012-11-16 09:02:50 +0100153 // test basic computations
154 {
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100155 DenseMatrix refM1 = DenseMatrix::Zero(rows, cols);
156 DenseMatrix refM2 = DenseMatrix::Zero(rows, cols);
157 DenseMatrix refM3 = DenseMatrix::Zero(rows, cols);
158 DenseMatrix refM4 = DenseMatrix::Zero(rows, cols);
159 SparseMatrixType m1(rows, cols);
160 SparseMatrixType m2(rows, cols);
161 SparseMatrixType m3(rows, cols);
162 SparseMatrixType m4(rows, cols);
Gael Guennebaud4e602832012-11-16 09:02:50 +0100163 initSparse<Scalar>(density, refM1, m1);
164 initSparse<Scalar>(density, refM2, m2);
165 initSparse<Scalar>(density, refM3, m3);
166 initSparse<Scalar>(density, refM4, m4);
167
Gael Guennebaud2c1b56f2016-05-31 10:56:53 +0200168 if(internal::random<bool>())
169 m1.makeCompressed();
170
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100171 Index m1_nnz = m1.nonZeros();
172
Gael Guennebaud4aac8722014-07-22 12:54:03 +0200173 VERIFY_IS_APPROX(m1*s1, refM1*s1);
Gael Guennebaud4e602832012-11-16 09:02:50 +0100174 VERIFY_IS_APPROX(m1+m2, refM1+refM2);
175 VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
176 VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
177 VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100178 VERIFY_IS_APPROX(m4=m1/s1, refM1/s1);
179 VERIFY_IS_EQUAL(m4.nonZeros(), m1_nnz);
Gael Guennebaud4e602832012-11-16 09:02:50 +0100180
Gael Guennebaud4e602832012-11-16 09:02:50 +0100181 if(SparseMatrixType::IsRowMajor)
182 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
183 else
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100184 VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.col(0)), refM1.col(0).dot(refM2.col(0)));
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100185
Gael Guennebaud3573a102014-02-17 13:46:17 +0100186 DenseVector rv = DenseVector::Random(m1.cols());
187 DenseVector cv = DenseVector::Random(m1.rows());
188 Index r = internal::random<Index>(0,m1.rows()-2);
189 Index c = internal::random<Index>(0,m1.cols()-1);
190 VERIFY_IS_APPROX(( m1.template block<1,Dynamic>(r,0,1,m1.cols()).dot(rv)) , refM1.row(r).dot(rv));
191 VERIFY_IS_APPROX(m1.row(r).dot(rv), refM1.row(r).dot(rv));
192 VERIFY_IS_APPROX(m1.col(c).dot(cv), refM1.col(c).dot(cv));
Gael Guennebaud4e602832012-11-16 09:02:50 +0100193
194 VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
195 VERIFY_IS_APPROX(m1.real(), refM1.real());
196
197 refM4.setRandom();
198 // sparse cwise* dense
199 VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
Gael Guennebaud90275082015-11-04 17:42:07 +0100200 // dense cwise* sparse
201 VERIFY_IS_APPROX(refM4.cwiseProduct(m3), refM4.cwiseProduct(refM3));
Gael Guennebaud4e602832012-11-16 09:02:50 +0100202// VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
203
Gael Guennebaudeec0dfd2018-10-10 22:50:15 +0200204 // mixed sparse-dense
Gael Guennebaud15084cf2016-01-29 22:09:45 +0100205 VERIFY_IS_APPROX(refM4 + m3, refM4 + refM3);
206 VERIFY_IS_APPROX(m3 + refM4, refM3 + refM4);
207 VERIFY_IS_APPROX(refM4 - m3, refM4 - refM3);
208 VERIFY_IS_APPROX(m3 - refM4, refM3 - refM4);
Gael Guennebaud71362672016-12-27 16:34:30 +0100209 VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + RealScalar(0.5)*m3).eval(), RealScalar(0.5)*refM4 + RealScalar(0.5)*refM3);
210 VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + m3*RealScalar(0.5)).eval(), RealScalar(0.5)*refM4 + RealScalar(0.5)*refM3);
211 VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + m3.cwiseProduct(m3)).eval(), RealScalar(0.5)*refM4 + refM3.cwiseProduct(refM3));
212
213 VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + RealScalar(0.5)*m3).eval(), RealScalar(0.5)*refM4 + RealScalar(0.5)*refM3);
214 VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + m3*RealScalar(0.5)).eval(), RealScalar(0.5)*refM4 + RealScalar(0.5)*refM3);
215 VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + (m3+m3)).eval(), RealScalar(0.5)*refM4 + (refM3+refM3));
216 VERIFY_IS_APPROX(((refM3+m3)+RealScalar(0.5)*m3).eval(), RealScalar(0.5)*refM3 + (refM3+refM3));
217 VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + (refM3+m3)).eval(), RealScalar(0.5)*refM4 + (refM3+refM3));
218 VERIFY_IS_APPROX((RealScalar(0.5)*refM4 + (m3+refM3)).eval(), RealScalar(0.5)*refM4 + (refM3+refM3));
219
Gael Guennebaud15084cf2016-01-29 22:09:45 +0100220
Gael Guennebaud2c1b56f2016-05-31 10:56:53 +0200221 VERIFY_IS_APPROX(m1.sum(), refM1.sum());
222
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100223 m4 = m1; refM4 = m4;
224
Gael Guennebaud2c1b56f2016-05-31 10:56:53 +0200225 VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100226 VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
Gael Guennebaud2c1b56f2016-05-31 10:56:53 +0200227 VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100228 VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
Gael Guennebaud2c1b56f2016-05-31 10:56:53 +0200229
230 VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
231 VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
232
Gael Guennebaudeec0dfd2018-10-10 22:50:15 +0200233 refM3 = refM1;
234
235 VERIFY_IS_APPROX(refM1+=m2, refM3+=refM2);
236 VERIFY_IS_APPROX(refM1-=m2, refM3-=refM2);
237
238 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1 =m2+refM4, refM3 =refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,10);
239 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1+=m2+refM4, refM3+=refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
240 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1-=m2+refM4, refM3-=refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
241 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1 =refM4+m2, refM3 =refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
242 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1+=refM4+m2, refM3+=refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
243 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1-=refM4+m2, refM3-=refM2+refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
244
245 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1 =m2-refM4, refM3 =refM2-refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,20);
246 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1+=m2-refM4, refM3+=refM2-refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
247 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1-=m2-refM4, refM3-=refM2-refM4); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
248 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1 =refM4-m2, refM3 =refM4-refM2); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
249 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1+=refM4-m2, refM3+=refM4-refM2); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
250 g_dense_op_sparse_count=0; VERIFY_IS_APPROX(refM1-=refM4-m2, refM3-=refM4-refM2); VERIFY_IS_EQUAL(g_dense_op_sparse_count,1);
251 refM3 = m3;
252
Gael Guennebaudba3f9772017-01-23 22:06:08 +0100253 if (rows>=2 && cols>=2)
254 {
255 VERIFY_RAISES_ASSERT( m1 += m1.innerVector(0) );
256 VERIFY_RAISES_ASSERT( m1 -= m1.innerVector(0) );
257 VERIFY_RAISES_ASSERT( refM1 -= m1.innerVector(0) );
258 VERIFY_RAISES_ASSERT( refM1 += m1.innerVector(0) );
259 }
Gael Guennebaud26a2c6f2017-12-14 15:11:04 +0100260 m1 = m4; refM1 = refM4;
Gael Guennebaudba3f9772017-01-23 22:06:08 +0100261
Gael Guennebaud4e602832012-11-16 09:02:50 +0100262 // test aliasing
263 VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100264 VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
265 m1 = m4; refM1 = refM4;
Gael Guennebaud4e602832012-11-16 09:02:50 +0100266 VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100267 VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
268 m1 = m4; refM1 = refM4;
Gael Guennebaud4e602832012-11-16 09:02:50 +0100269 VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100270 VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
271 m1 = m4; refM1 = refM4;
Gael Guennebaud4e602832012-11-16 09:02:50 +0100272 VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
Gael Guennebaudc86911a2017-01-30 13:38:24 +0100273 VERIFY_IS_EQUAL(m1.nonZeros(), m1_nnz);
274 m1 = m4; refM1 = refM4;
Gael Guennebaud7e029d12016-08-29 12:06:37 +0200275
276 if(m1.isCompressed())
277 {
278 VERIFY_IS_APPROX(m1.coeffs().sum(), m1.sum());
279 m1.coeffs() += s1;
280 for(Index j = 0; j<m1.outerSize(); ++j)
281 for(typename SparseMatrixType::InnerIterator it(m1,j); it; ++it)
282 refM1(it.row(), it.col()) += s1;
283 VERIFY_IS_APPROX(m1, refM1);
284 }
Gael Guennebaud2e334f52016-11-14 18:47:02 +0100285
286 // and/or
287 {
288 typedef SparseMatrix<bool, SparseMatrixType::Options, typename SparseMatrixType::StorageIndex> SpBool;
289 SpBool mb1 = m1.real().template cast<bool>();
290 SpBool mb2 = m2.real().template cast<bool>();
291 VERIFY_IS_EQUAL(mb1.template cast<int>().sum(), refM1.real().template cast<bool>().count());
292 VERIFY_IS_EQUAL((mb1 && mb2).template cast<int>().sum(), (refM1.real().template cast<bool>() && refM2.real().template cast<bool>()).count());
293 VERIFY_IS_EQUAL((mb1 || mb2).template cast<int>().sum(), (refM1.real().template cast<bool>() || refM2.real().template cast<bool>()).count());
294 SpBool mb3 = mb1 && mb2;
295 if(mb1.coeffs().all() && mb2.coeffs().all())
296 {
297 VERIFY_IS_EQUAL(mb3.nonZeros(), (refM1.real().template cast<bool>() && refM2.real().template cast<bool>()).count());
298 }
299 }
Gael Guennebaud4e602832012-11-16 09:02:50 +0100300 }
301
Gael Guennebaudeedb87f2016-11-14 14:05:53 +0100302 // test reverse iterators
303 {
304 DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
305 SparseMatrixType m2(rows, cols);
306 initSparse<Scalar>(density, refMat2, m2);
307 std::vector<Scalar> ref_value(m2.innerSize());
308 std::vector<Index> ref_index(m2.innerSize());
309 if(internal::random<bool>())
310 m2.makeCompressed();
311 for(Index j = 0; j<m2.outerSize(); ++j)
312 {
313 Index count_forward = 0;
314
315 for(typename SparseMatrixType::InnerIterator it(m2,j); it; ++it)
316 {
317 ref_value[ref_value.size()-1-count_forward] = it.value();
318 ref_index[ref_index.size()-1-count_forward] = it.index();
319 count_forward++;
320 }
321 Index count_reverse = 0;
322 for(typename SparseMatrixType::ReverseInnerIterator it(m2,j); it; --it)
323 {
324 VERIFY_IS_APPROX( std::abs(ref_value[ref_value.size()-count_forward+count_reverse])+1, std::abs(it.value())+1);
325 VERIFY_IS_EQUAL( ref_index[ref_index.size()-count_forward+count_reverse] , it.index());
326 count_reverse++;
327 }
328 VERIFY_IS_EQUAL(count_forward, count_reverse);
329 }
330 }
331
Gael Guennebaud4e602832012-11-16 09:02:50 +0100332 // test transpose
333 {
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100334 DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
335 SparseMatrixType m2(rows, cols);
Gael Guennebaud4e602832012-11-16 09:02:50 +0100336 initSparse<Scalar>(density, refMat2, m2);
337 VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
338 VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
339
340 VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
Gael Guennebaudff46ec02014-09-22 23:33:28 +0200341
342 // check isApprox handles opposite storage order
343 typename Transpose<SparseMatrixType>::PlainObject m3(m2);
344 VERIFY(m2.isApprox(m3));
Gael Guennebaud4e602832012-11-16 09:02:50 +0100345 }
346
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000347 // test prune
348 {
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100349 SparseMatrixType m2(rows, cols);
350 DenseMatrix refM2(rows, cols);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000351 refM2.setZero();
352 int countFalseNonZero = 0;
353 int countTrueNonZero = 0;
Gael Guennebauda44d91a2015-10-13 10:53:38 +0200354 m2.reserve(VectorXi::Constant(m2.outerSize(), int(m2.innerSize())));
355 for (Index j=0; j<m2.cols(); ++j)
Gael Guennebaud28293142009-05-04 14:25:12 +0000356 {
Gael Guennebauda44d91a2015-10-13 10:53:38 +0200357 for (Index i=0; i<m2.rows(); ++i)
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000358 {
Benoit Jacob47160402010-10-25 10:15:22 -0400359 float x = internal::random<float>(0,1);
Christoph Hertzbergdacb4692016-05-05 13:35:45 +0200360 if (x<0.1f)
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000361 {
362 // do nothing
363 }
Christoph Hertzbergdacb4692016-05-05 13:35:45 +0200364 else if (x<0.5f)
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000365 {
366 countFalseNonZero++;
Gael Guennebauda44d91a2015-10-13 10:53:38 +0200367 m2.insert(i,j) = Scalar(0);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000368 }
369 else
370 {
371 countTrueNonZero++;
Gael Guennebauda44d91a2015-10-13 10:53:38 +0200372 m2.insert(i,j) = Scalar(1);
373 refM2(i,j) = Scalar(1);
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000374 }
375 }
Gael Guennebaud28293142009-05-04 14:25:12 +0000376 }
Gael Guennebauda44d91a2015-10-13 10:53:38 +0200377 if(internal::random<bool>())
378 m2.makeCompressed();
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000379 VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
Gael Guennebauda44d91a2015-10-13 10:53:38 +0200380 if(countTrueNonZero>0)
381 VERIFY_IS_APPROX(m2, refM2);
Hauke Heibeld204ec42010-11-02 14:33:33 +0100382 m2.prune(Scalar(1));
Gael Guennebaud52cf07d2009-01-21 18:46:04 +0000383 VERIFY(countTrueNonZero==m2.nonZeros());
384 VERIFY_IS_APPROX(m2, refM2);
385 }
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100386
Gael Guennebaud87138072012-01-28 11:13:59 +0100387 // test setFromTriplets
388 {
Christoph Hertzberge8cdbed2014-12-04 22:48:53 +0100389 typedef Triplet<Scalar,StorageIndex> TripletType;
Gael Guennebaud87138072012-01-28 11:13:59 +0100390 std::vector<TripletType> triplets;
Gael Guennebaudb47ef142014-07-08 16:47:11 +0200391 Index ntriplets = rows*cols;
Gael Guennebaud87138072012-01-28 11:13:59 +0100392 triplets.reserve(ntriplets);
Gael Guennebaudb4c79ee2015-10-13 11:30:41 +0200393 DenseMatrix refMat_sum = DenseMatrix::Zero(rows,cols);
394 DenseMatrix refMat_prod = DenseMatrix::Zero(rows,cols);
395 DenseMatrix refMat_last = DenseMatrix::Zero(rows,cols);
396
Gael Guennebaudb47ef142014-07-08 16:47:11 +0200397 for(Index i=0;i<ntriplets;++i)
Gael Guennebaud87138072012-01-28 11:13:59 +0100398 {
Gael Guennebaudaa6c5162015-02-16 13:19:05 +0100399 StorageIndex r = internal::random<StorageIndex>(0,StorageIndex(rows-1));
400 StorageIndex c = internal::random<StorageIndex>(0,StorageIndex(cols-1));
Gael Guennebaud87138072012-01-28 11:13:59 +0100401 Scalar v = internal::random<Scalar>();
Gael Guennebaud18e3ac02012-01-31 09:14:01 +0100402 triplets.push_back(TripletType(r,c,v));
Gael Guennebaudb4c79ee2015-10-13 11:30:41 +0200403 refMat_sum(r,c) += v;
404 if(std::abs(refMat_prod(r,c))==0)
405 refMat_prod(r,c) = v;
406 else
407 refMat_prod(r,c) *= v;
408 refMat_last(r,c) = v;
Gael Guennebaud87138072012-01-28 11:13:59 +0100409 }
410 SparseMatrixType m(rows,cols);
411 m.setFromTriplets(triplets.begin(), triplets.end());
Gael Guennebaudb4c79ee2015-10-13 11:30:41 +0200412 VERIFY_IS_APPROX(m, refMat_sum);
413
414 m.setFromTriplets(triplets.begin(), triplets.end(), std::multiplies<Scalar>());
415 VERIFY_IS_APPROX(m, refMat_prod);
416#if (defined(__cplusplus) && __cplusplus >= 201103L)
417 m.setFromTriplets(triplets.begin(), triplets.end(), [] (Scalar,Scalar b) { return b; });
418 VERIFY_IS_APPROX(m, refMat_last);
419#endif
Gael Guennebaud87138072012-01-28 11:13:59 +0100420 }
Gael Guennebaud3af29ca2015-02-09 10:23:45 +0100421
422 // test Map
423 {
424 DenseMatrix refMat2(rows, cols), refMat3(rows, cols);
425 SparseMatrixType m2(rows, cols), m3(rows, cols);
426 initSparse<Scalar>(density, refMat2, m2);
427 initSparse<Scalar>(density, refMat3, m3);
428 {
429 Map<SparseMatrixType> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr());
430 Map<SparseMatrixType> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr());
431 VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
432 VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
433 }
434 {
Gael Guennebaudfe513192015-02-13 10:03:53 +0100435 MappedSparseMatrix<Scalar,SparseMatrixType::Options,StorageIndex> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr());
436 MappedSparseMatrix<Scalar,SparseMatrixType::Options,StorageIndex> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr());
Gael Guennebaud3af29ca2015-02-09 10:23:45 +0100437 VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
438 VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3);
439 }
Gael Guennebaud757971e2016-07-26 09:40:19 +0200440
441 Index i = internal::random<Index>(0,rows-1);
442 Index j = internal::random<Index>(0,cols-1);
443 m2.coeffRef(i,j) = 123;
444 if(internal::random<bool>())
445 m2.makeCompressed();
446 Map<SparseMatrixType> mapMat2(rows, cols, m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr());
447 VERIFY_IS_EQUAL(m2.coeff(i,j),Scalar(123));
448 VERIFY_IS_EQUAL(mapMat2.coeff(i,j),Scalar(123));
449 mapMat2.coeffRef(i,j) = -123;
450 VERIFY_IS_EQUAL(m2.coeff(i,j),Scalar(-123));
Gael Guennebaud3af29ca2015-02-09 10:23:45 +0100451 }
Gael Guennebaud87138072012-01-28 11:13:59 +0100452
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100453 // test triangularView
454 {
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100455 DenseMatrix refMat2(rows, cols), refMat3(rows, cols);
456 SparseMatrixType m2(rows, cols), m3(rows, cols);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100457 initSparse<Scalar>(density, refMat2, m2);
458 refMat3 = refMat2.template triangularView<Lower>();
459 m3 = m2.template triangularView<Lower>();
460 VERIFY_IS_APPROX(m3, refMat3);
461
462 refMat3 = refMat2.template triangularView<Upper>();
463 m3 = m2.template triangularView<Upper>();
464 VERIFY_IS_APPROX(m3, refMat3);
465
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100466 {
467 refMat3 = refMat2.template triangularView<UnitUpper>();
468 m3 = m2.template triangularView<UnitUpper>();
469 VERIFY_IS_APPROX(m3, refMat3);
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100470
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100471 refMat3 = refMat2.template triangularView<UnitLower>();
472 m3 = m2.template triangularView<UnitLower>();
473 VERIFY_IS_APPROX(m3, refMat3);
474 }
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200475
476 refMat3 = refMat2.template triangularView<StrictlyUpper>();
477 m3 = m2.template triangularView<StrictlyUpper>();
478 VERIFY_IS_APPROX(m3, refMat3);
479
480 refMat3 = refMat2.template triangularView<StrictlyLower>();
481 m3 = m2.template triangularView<StrictlyLower>();
482 VERIFY_IS_APPROX(m3, refMat3);
Gael Guennebaudf6b1dee2015-11-04 17:02:32 +0100483
Gael Guennebaudba3f9772017-01-23 22:06:08 +0100484 // check sparse-triangular to dense
Gael Guennebaudf6b1dee2015-11-04 17:02:32 +0100485 refMat3 = m2.template triangularView<StrictlyUpper>();
486 VERIFY_IS_APPROX(refMat3, DenseMatrix(refMat2.template triangularView<StrictlyUpper>()));
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100487 }
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200488
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100489 // test selfadjointView
Gael Guennebaud9353bba2011-12-04 14:39:24 +0100490 if(!SparseMatrixType::IsRowMajor)
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100491 {
492 DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
493 SparseMatrixType m2(rows, rows), m3(rows, rows);
494 initSparse<Scalar>(density, refMat2, m2);
495 refMat3 = refMat2.template selfadjointView<Lower>();
496 m3 = m2.template selfadjointView<Lower>();
497 VERIFY_IS_APPROX(m3, refMat3);
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100498
Gael Guennebaud11b492e2016-12-14 17:53:47 +0100499 refMat3 += refMat2.template selfadjointView<Lower>();
500 m3 += m2.template selfadjointView<Lower>();
501 VERIFY_IS_APPROX(m3, refMat3);
502
503 refMat3 -= refMat2.template selfadjointView<Lower>();
504 m3 -= m2.template selfadjointView<Lower>();
505 VERIFY_IS_APPROX(m3, refMat3);
506
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100507 // selfadjointView only works for square matrices:
508 SparseMatrixType m4(rows, rows+1);
509 VERIFY_RAISES_ASSERT(m4.template selfadjointView<Lower>());
510 VERIFY_RAISES_ASSERT(m4.template selfadjointView<Upper>());
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100511 }
512
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200513 // test sparseView
514 {
515 DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
516 SparseMatrixType m2(rows, rows);
Gael Guennebaud9a3ec632010-11-15 14:14:05 +0100517 initSparse<Scalar>(density, refMat2, m2);
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200518 VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
Gael Guennebaud8456bbb2016-05-18 16:53:28 +0200519
520 // sparse view on expressions:
521 VERIFY_IS_APPROX((s1*m2).eval(), (s1*refMat2).sparseView().eval());
522 VERIFY_IS_APPROX((m2+m2).eval(), (refMat2+refMat2).sparseView().eval());
523 VERIFY_IS_APPROX((m2*m2).eval(), (refMat2.lazyProduct(refMat2)).sparseView().eval());
524 VERIFY_IS_APPROX((m2*m2).eval(), (refMat2*refMat2).sparseView().eval());
Gael Guennebaudfa6d36e2010-07-22 15:57:01 +0200525 }
Gael Guennebaud82f9aa12011-12-04 21:49:21 +0100526
527 // test diagonal
528 {
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100529 DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
530 SparseMatrixType m2(rows, cols);
Gael Guennebaud82f9aa12011-12-04 21:49:21 +0100531 initSparse<Scalar>(density, refMat2, m2);
532 VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
Gael Guennebaud296d24b2017-01-25 17:39:01 +0100533 DenseVector d = m2.diagonal();
534 VERIFY_IS_APPROX(d, refMat2.diagonal().eval());
535 d = m2.diagonal().array();
536 VERIFY_IS_APPROX(d, refMat2.diagonal().eval());
Gael Guennebaudb26e6972014-12-01 14:41:39 +0100537 VERIFY_IS_APPROX(const_cast<const SparseMatrixType&>(m2).diagonal(), refMat2.diagonal().eval());
538
539 initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag);
540 m2.diagonal() += refMat2.diagonal();
541 refMat2.diagonal() += refMat2.diagonal();
542 VERIFY_IS_APPROX(m2, refMat2);
Gael Guennebaud82f9aa12011-12-04 21:49:21 +0100543 }
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200544
Gael Guennebaud62f21e22015-06-24 17:55:00 +0200545 // test diagonal to sparse
546 {
547 DenseVector d = DenseVector::Random(rows);
548 DenseMatrix refMat2 = d.asDiagonal();
Gael Guennebaudf489f442019-01-28 17:29:50 +0100549 SparseMatrixType m2;
Gael Guennebaud62f21e22015-06-24 17:55:00 +0200550 m2 = d.asDiagonal();
551 VERIFY_IS_APPROX(m2, refMat2);
Gael Guennebaud4c8cd132015-06-24 18:11:06 +0200552 SparseMatrixType m3(d.asDiagonal());
553 VERIFY_IS_APPROX(m3, refMat2);
Gael Guennebaud62f21e22015-06-24 17:55:00 +0200554 refMat2 += d.asDiagonal();
555 m2 += d.asDiagonal();
556 VERIFY_IS_APPROX(m2, refMat2);
Gael Guennebaudf489f442019-01-28 17:29:50 +0100557 m2.setZero(); m2 += d.asDiagonal();
558 refMat2.setZero(); refMat2 += d.asDiagonal();
559 VERIFY_IS_APPROX(m2, refMat2);
560 m2.setZero(); m2 -= d.asDiagonal();
561 refMat2.setZero(); refMat2 -= d.asDiagonal();
562 VERIFY_IS_APPROX(m2, refMat2);
563
564 initSparse<Scalar>(density, refMat2, m2);
565 m2.makeCompressed();
566 m2 += d.asDiagonal();
567 refMat2 += d.asDiagonal();
568 VERIFY_IS_APPROX(m2, refMat2);
569
570 initSparse<Scalar>(density, refMat2, m2);
571 m2.makeCompressed();
572 VectorXi res(rows);
573 for(Index i=0; i<rows; ++i)
574 res(i) = internal::random<int>(0,3);
575 m2.reserve(res);
576 m2 -= d.asDiagonal();
577 refMat2 -= d.asDiagonal();
578 VERIFY_IS_APPROX(m2, refMat2);
Gael Guennebaud62f21e22015-06-24 17:55:00 +0200579 }
580
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200581 // test conservative resize
582 {
Christoph Hertzberge8cdbed2014-12-04 22:48:53 +0100583 std::vector< std::pair<StorageIndex,StorageIndex> > inc;
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100584 if(rows > 3 && cols > 2)
Christoph Hertzberge8cdbed2014-12-04 22:48:53 +0100585 inc.push_back(std::pair<StorageIndex,StorageIndex>(-3,-2));
586 inc.push_back(std::pair<StorageIndex,StorageIndex>(0,0));
587 inc.push_back(std::pair<StorageIndex,StorageIndex>(3,2));
588 inc.push_back(std::pair<StorageIndex,StorageIndex>(3,0));
589 inc.push_back(std::pair<StorageIndex,StorageIndex>(0,3));
Adam Shapiro2ac0b782021-02-23 21:32:39 +0000590 inc.push_back(std::pair<StorageIndex,StorageIndex>(0,-1));
591 inc.push_back(std::pair<StorageIndex,StorageIndex>(-1,0));
592 inc.push_back(std::pair<StorageIndex,StorageIndex>(-1,-1));
593
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200594 for(size_t i = 0; i< inc.size(); i++) {
Christoph Hertzberge8cdbed2014-12-04 22:48:53 +0100595 StorageIndex incRows = inc[i].first;
596 StorageIndex incCols = inc[i].second;
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200597 SparseMatrixType m1(rows, cols);
598 DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
599 initSparse<Scalar>(density, refMat1, m1);
Adam Shapiro2ac0b782021-02-23 21:32:39 +0000600
601 SparseMatrixType m2 = m1;
602 m2.makeCompressed();
603
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200604 m1.conservativeResize(rows+incRows, cols+incCols);
Adam Shapiro2ac0b782021-02-23 21:32:39 +0000605 m2.conservativeResize(rows+incRows, cols+incCols);
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200606 refMat1.conservativeResize(rows+incRows, cols+incCols);
607 if (incRows > 0) refMat1.bottomRows(incRows).setZero();
608 if (incCols > 0) refMat1.rightCols(incCols).setZero();
Adam Shapiro2ac0b782021-02-23 21:32:39 +0000609
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200610 VERIFY_IS_APPROX(m1, refMat1);
Adam Shapiro2ac0b782021-02-23 21:32:39 +0000611 VERIFY_IS_APPROX(m2, refMat1);
612
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200613 // Insert new values
614 if (incRows > 0)
Gael Guennebaud7ee378d2013-07-12 16:40:02 +0200615 m1.insert(m1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200616 if (incCols > 0)
Gael Guennebaud7ee378d2013-07-12 16:40:02 +0200617 m1.insert(0, m1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
Adam Shapiro2ac0b782021-02-23 21:32:39 +0000618
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200619 VERIFY_IS_APPROX(m1, refMat1);
Adam Shapiro2ac0b782021-02-23 21:32:39 +0000620
621
Benjamin Piwowarski6bf49ce2012-07-19 00:07:06 +0200622 }
623 }
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200624
625 // test Identity matrix
626 {
627 DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows);
628 SparseMatrixType m1(rows, rows);
629 m1.setIdentity();
630 VERIFY_IS_APPROX(m1, refMat1);
Gael Guennebaud8a211bb2015-10-25 22:01:58 +0100631 for(int k=0; k<rows*rows/4; ++k)
632 {
633 Index i = internal::random<Index>(0,rows-1);
634 Index j = internal::random<Index>(0,rows-1);
Gael Guennebaud73f692d2015-10-27 11:01:37 +0100635 Scalar v = internal::random<Scalar>();
Gael Guennebaud8a211bb2015-10-25 22:01:58 +0100636 m1.coeffRef(i,j) = v;
637 refMat1.coeffRef(i,j) = v;
638 VERIFY_IS_APPROX(m1, refMat1);
639 if(internal::random<Index>(0,10)<2)
640 m1.makeCompressed();
641 }
642 m1.setIdentity();
643 refMat1.setIdentity();
644 VERIFY_IS_APPROX(m1, refMat1);
Desire NUENTSA4cd82452013-06-11 14:42:29 +0200645 }
Gael Guennebaudec469702016-02-01 15:04:33 +0100646
647 // test array/vector of InnerIterator
648 {
649 typedef typename SparseMatrixType::InnerIterator IteratorType;
650
651 DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols);
652 SparseMatrixType m2(rows, cols);
653 initSparse<Scalar>(density, refMat2, m2);
654 IteratorType static_array[2];
655 static_array[0] = IteratorType(m2,0);
656 static_array[1] = IteratorType(m2,m2.outerSize()-1);
657 VERIFY( static_array[0] || m2.innerVector(static_array[0].outer()).nonZeros() == 0 );
658 VERIFY( static_array[1] || m2.innerVector(static_array[1].outer()).nonZeros() == 0 );
659 if(static_array[0] && static_array[1])
660 {
661 ++(static_array[1]);
662 static_array[1] = IteratorType(m2,0);
663 VERIFY( static_array[1] );
664 VERIFY( static_array[1].index() == static_array[0].index() );
665 VERIFY( static_array[1].outer() == static_array[0].outer() );
666 VERIFY( static_array[1].value() == static_array[0].value() );
667 }
668
669 std::vector<IteratorType> iters(2);
670 iters[0] = IteratorType(m2,0);
671 iters[1] = IteratorType(m2,m2.outerSize()-1);
672 }
Gael Guennebaud25424d92020-08-26 12:32:20 +0200673
674 // test reserve with empty rows/columns
675 {
676 SparseMatrixType m1(0,cols);
677 m1.reserve(ArrayXi::Constant(m1.outerSize(),1));
678 SparseMatrixType m2(rows,0);
679 m2.reserve(ArrayXi::Constant(m2.outerSize(),1));
680 }
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000681}
682
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100683
Christoph Hertzbergc5a37772014-10-31 17:19:05 +0100684template<typename SparseMatrixType>
Christoph Hertzberge8cdbed2014-12-04 22:48:53 +0100685void big_sparse_triplet(Index rows, Index cols, double density) {
686 typedef typename SparseMatrixType::StorageIndex StorageIndex;
687 typedef typename SparseMatrixType::Scalar Scalar;
688 typedef Triplet<Scalar,Index> TripletType;
689 std::vector<TripletType> triplets;
690 double nelements = density * rows*cols;
691 VERIFY(nelements>=0 && nelements < NumTraits<StorageIndex>::highest());
692 Index ntriplets = Index(nelements);
693 triplets.reserve(ntriplets);
694 Scalar sum = Scalar(0);
695 for(Index i=0;i<ntriplets;++i)
696 {
697 Index r = internal::random<Index>(0,rows-1);
698 Index c = internal::random<Index>(0,cols-1);
Gael Guennebaudcd25b532018-12-08 00:13:37 +0100699 // use positive values to prevent numerical cancellation errors in sum
700 Scalar v = numext::abs(internal::random<Scalar>());
Christoph Hertzberge8cdbed2014-12-04 22:48:53 +0100701 triplets.push_back(TripletType(r,c,v));
702 sum += v;
703 }
704 SparseMatrixType m(rows,cols);
705 m.setFromTriplets(triplets.begin(), triplets.end());
706 VERIFY(m.nonZeros() <= ntriplets);
707 VERIFY_IS_APPROX(sum, m.sum());
Christoph Hertzbergc5a37772014-10-31 17:19:05 +0100708}
709
Gael Guennebauddff3a922018-07-17 15:52:58 +0200710template<int>
711void bug1105()
712{
713 // Regression test for bug 1105
714 int n = Eigen::internal::random<int>(200,600);
715 SparseMatrix<std::complex<double>,0, long> mat(n, n);
716 std::complex<double> val;
717
718 for(int i=0; i<n; ++i)
719 {
720 mat.coeffRef(i, i%(n/10)) = val;
721 VERIFY(mat.data().allocatedSize()<20*n);
722 }
723}
Christoph Hertzbergc5a37772014-10-31 17:19:05 +0100724
Gael Guennebaud8214cf12018-10-11 10:27:23 +0200725#ifndef EIGEN_SPARSE_TEST_INCLUDED_FROM_SPARSE_EXTRA
726
Gael Guennebaud82f0ce22018-07-17 14:46:15 +0200727EIGEN_DECLARE_TEST(sparse_basic)
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000728{
Rasmus Munk Larsen954b4ca2018-10-22 13:48:56 -0700729 g_dense_op_sparse_count = 0; // Suppresses compiler warning.
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000730 for(int i = 0; i < g_repeat; i++) {
Gael Guennebaudc43154b2015-03-04 10:16:46 +0100731 int r = Eigen::internal::random<int>(1,200), c = Eigen::internal::random<int>(1,200);
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100732 if(Eigen::internal::random<int>(0,4) == 0) {
733 r = c; // check square matrices in 25% of tries
734 }
735 EIGEN_UNUSED_VARIABLE(r+c);
736 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(1, 1)) ));
Gael Guennebaud91fe1502011-06-07 11:28:16 +0200737 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
Christoph Hertzberg0833b822014-10-31 17:12:13 +0100738 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(r, c)) ));
739 CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(r, c)) ));
740 CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(r, c)) ));
Gael Guennebaudd7698c12015-03-19 15:11:05 +0100741 CALL_SUBTEST_5(( sparse_basic(SparseMatrix<double,ColMajor,long int>(r, c)) ));
742 CALL_SUBTEST_5(( sparse_basic(SparseMatrix<double,RowMajor,long int>(r, c)) ));
Gael Guennebaud47e89022013-06-10 10:34:03 +0200743
Gael Guennebaudc43154b2015-03-04 10:16:46 +0100744 r = Eigen::internal::random<int>(1,100);
745 c = Eigen::internal::random<int>(1,100);
746 if(Eigen::internal::random<int>(0,4) == 0) {
747 r = c; // check square matrices in 25% of tries
748 }
749
Gael Guennebaudd7698c12015-03-19 15:11:05 +0100750 CALL_SUBTEST_6(( sparse_basic(SparseMatrix<double,ColMajor,short int>(short(r), short(c))) ));
751 CALL_SUBTEST_6(( sparse_basic(SparseMatrix<double,RowMajor,short int>(short(r), short(c))) ));
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000752 }
Christoph Hertzbergc5a37772014-10-31 17:19:05 +0100753
754 // Regression test for bug 900: (manually insert higher values here, if you have enough RAM):
755 CALL_SUBTEST_3((big_sparse_triplet<SparseMatrix<float, RowMajor, int> >(10000, 10000, 0.125)));
756 CALL_SUBTEST_4((big_sparse_triplet<SparseMatrix<double, ColMajor, long int> >(10000, 10000, 0.125)));
Gael Guennebaudbfd6ee62015-11-06 15:05:37 +0100757
Gael Guennebauddff3a922018-07-17 15:52:58 +0200758 CALL_SUBTEST_7( bug1105<0>() );
Gael Guennebaud86ccd992008-11-05 13:47:55 +0000759}
Gael Guennebaud8214cf12018-10-11 10:27:23 +0200760#endif