Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 1 | // This file is part of Eigen, a lightweight C++ template library |
Benoit Jacob | 6347b1d | 2009-05-22 20:25:33 +0200 | [diff] [blame] | 2 | // for linear algebra. |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 3 | // |
Gael Guennebaud | 22c7609 | 2011-03-22 11:58:22 +0100 | [diff] [blame] | 4 | // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr> |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 5 | // Copyright (C) 2008 Daniel Gomez Ferro <dgomezferro@gmail.com> |
Desire NUENTSA | 4cd8245 | 2013-06-11 14:42:29 +0200 | [diff] [blame] | 6 | // Copyright (C) 2013 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr> |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 7 | // |
Benoit Jacob | 69124cf | 2012-07-13 14:42:47 -0400 | [diff] [blame] | 8 | // 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 Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 11 | |
Gael Guennebaud | c43154b | 2015-03-04 10:16:46 +0100 | [diff] [blame] | 12 | static long g_realloc_count = 0; |
| 13 | #define EIGEN_SPARSE_COMPRESSED_STORAGE_REALLOCATE_PLUGIN g_realloc_count++; |
| 14 | |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 15 | #include "sparse.h" |
| 16 | |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 17 | template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref) |
| 18 | { |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 19 | typedef typename SparseMatrixType::StorageIndex StorageIndex; |
| 20 | typedef Matrix<StorageIndex,2,1> Vector2; |
Gael Guennebaud | 6d1f5db | 2013-07-10 23:48:26 +0200 | [diff] [blame] | 21 | |
Gael Guennebaud | fc202ba | 2015-02-13 18:57:41 +0100 | [diff] [blame] | 22 | const Index rows = ref.rows(); |
| 23 | const Index cols = ref.cols(); |
Gael Guennebaud | 0a537cb | 2016-02-12 15:58:31 +0100 | [diff] [blame] | 24 | //const Index inner = ref.innerSize(); |
| 25 | //const Index outer = ref.outerSize(); |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 26 | |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 27 | typedef typename SparseMatrixType::Scalar Scalar; |
| 28 | enum { Flags = SparseMatrixType::Flags }; |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 29 | |
Gael Guennebaud | 42e2578 | 2011-08-19 14:18:05 +0200 | [diff] [blame] | 30 | double density = (std::max)(8./(rows*cols), 0.01); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 31 | typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix; |
| 32 | typedef Matrix<Scalar,Dynamic,1> DenseVector; |
| 33 | Scalar eps = 1e-6; |
| 34 | |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 35 | Scalar s1 = internal::random<Scalar>(); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 36 | { |
Gael Guennebaud | d1d7a1a | 2013-06-23 19:11:32 +0200 | [diff] [blame] | 37 | SparseMatrixType m(rows, cols); |
| 38 | DenseMatrix refMat = DenseMatrix::Zero(rows, cols); |
| 39 | DenseVector vec1 = DenseVector::Random(rows); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 40 | |
Gael Guennebaud | 6d1f5db | 2013-07-10 23:48:26 +0200 | [diff] [blame] | 41 | std::vector<Vector2> zeroCoords; |
| 42 | std::vector<Vector2> nonzeroCoords; |
Gael Guennebaud | d1d7a1a | 2013-06-23 19:11:32 +0200 | [diff] [blame] | 43 | initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 44 | |
Gael Guennebaud | d1d7a1a | 2013-06-23 19:11:32 +0200 | [diff] [blame] | 45 | // test coeff and coeffRef |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 46 | for (std::size_t i=0; i<zeroCoords.size(); ++i) |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 47 | { |
Gael Guennebaud | d1d7a1a | 2013-06-23 19:11:32 +0200 | [diff] [blame] | 48 | VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps ); |
| 49 | if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value) |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 50 | VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[i].x(),zeroCoords[i].y()) = 5 ); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 51 | } |
Gael Guennebaud | d1d7a1a | 2013-06-23 19:11:32 +0200 | [diff] [blame] | 52 | VERIFY_IS_APPROX(m, refMat); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 53 | |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 54 | if(!nonzeroCoords.empty()) { |
| 55 | m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5); |
| 56 | refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5); |
| 57 | } |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 58 | |
Gael Guennebaud | d1d7a1a | 2013-06-23 19:11:32 +0200 | [diff] [blame] | 59 | VERIFY_IS_APPROX(m, refMat); |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 60 | |
Gael Guennebaud | a915f02 | 2013-06-28 16:16:02 +0200 | [diff] [blame] | 61 | // test assertion |
| 62 | VERIFY_RAISES_ASSERT( m.coeffRef(-1,1) = 0 ); |
| 63 | VERIFY_RAISES_ASSERT( m.coeffRef(0,m.cols()) = 0 ); |
Gael Guennebaud | d1d7a1a | 2013-06-23 19:11:32 +0200 | [diff] [blame] | 64 | } |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 65 | |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 66 | // test insert (inner random) |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 67 | { |
| 68 | DenseMatrix m1(rows,cols); |
| 69 | m1.setZero(); |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 70 | SparseMatrixType m2(rows,cols); |
Gael Guennebaud | c43154b | 2015-03-04 10:16:46 +0100 | [diff] [blame] | 71 | bool call_reserve = internal::random<int>()%2; |
| 72 | Index nnz = internal::random<int>(1,int(rows)/2); |
| 73 | if(call_reserve) |
| 74 | { |
| 75 | if(internal::random<int>()%2) |
| 76 | m2.reserve(VectorXi::Constant(m2.outerSize(), int(nnz))); |
| 77 | else |
| 78 | m2.reserve(m2.outerSize() * nnz); |
| 79 | } |
| 80 | g_realloc_count = 0; |
Gael Guennebaud | 6d1f5db | 2013-07-10 23:48:26 +0200 | [diff] [blame] | 81 | for (Index j=0; j<cols; ++j) |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 82 | { |
Gael Guennebaud | c43154b | 2015-03-04 10:16:46 +0100 | [diff] [blame] | 83 | for (Index k=0; k<nnz; ++k) |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 84 | { |
Gael Guennebaud | 6d1f5db | 2013-07-10 23:48:26 +0200 | [diff] [blame] | 85 | Index i = internal::random<Index>(0,rows-1); |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 86 | if (m1.coeff(i,j)==Scalar(0)) |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 87 | m2.insert(i,j) = m1(i,j) = internal::random<Scalar>(); |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 88 | } |
| 89 | } |
Gael Guennebaud | c43154b | 2015-03-04 10:16:46 +0100 | [diff] [blame] | 90 | |
| 91 | if(call_reserve && !SparseMatrixType::IsRowMajor) |
| 92 | { |
| 93 | VERIFY(g_realloc_count==0); |
| 94 | } |
| 95 | |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 96 | m2.finalize(); |
| 97 | VERIFY_IS_APPROX(m2,m1); |
| 98 | } |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 99 | |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 100 | // test insert (fully random) |
| 101 | { |
| 102 | DenseMatrix m1(rows,cols); |
| 103 | m1.setZero(); |
| 104 | SparseMatrixType m2(rows,cols); |
Gael Guennebaud | 4ca89f3 | 2011-12-02 19:00:16 +0100 | [diff] [blame] | 105 | if(internal::random<int>()%2) |
| 106 | m2.reserve(VectorXi::Constant(m2.outerSize(), 2)); |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 107 | for (int k=0; k<rows*cols; ++k) |
| 108 | { |
Gael Guennebaud | 6d1f5db | 2013-07-10 23:48:26 +0200 | [diff] [blame] | 109 | Index i = internal::random<Index>(0,rows-1); |
| 110 | Index j = internal::random<Index>(0,cols-1); |
Gael Guennebaud | 4ca89f3 | 2011-12-02 19:00:16 +0100 | [diff] [blame] | 111 | if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2)) |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 112 | m2.insert(i,j) = m1(i,j) = internal::random<Scalar>(); |
Gael Guennebaud | 4ca89f3 | 2011-12-02 19:00:16 +0100 | [diff] [blame] | 113 | else |
| 114 | { |
| 115 | Scalar v = internal::random<Scalar>(); |
| 116 | m2.coeffRef(i,j) += v; |
| 117 | m1(i,j) += v; |
| 118 | } |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 119 | } |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 120 | VERIFY_IS_APPROX(m2,m1); |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 121 | } |
Gael Guennebaud | 7706baf | 2011-09-08 13:42:54 +0200 | [diff] [blame] | 122 | |
| 123 | // test insert (un-compressed) |
| 124 | for(int mode=0;mode<4;++mode) |
| 125 | { |
| 126 | DenseMatrix m1(rows,cols); |
| 127 | m1.setZero(); |
| 128 | SparseMatrixType m2(rows,cols); |
Gael Guennebaud | b47ef14 | 2014-07-08 16:47:11 +0200 | [diff] [blame] | 129 | VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? int(m2.innerSize()) : std::max<int>(1,int(m2.innerSize())/8))); |
Gael Guennebaud | 7706baf | 2011-09-08 13:42:54 +0200 | [diff] [blame] | 130 | m2.reserve(r); |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 131 | for (Index k=0; k<rows*cols; ++k) |
Gael Guennebaud | 7706baf | 2011-09-08 13:42:54 +0200 | [diff] [blame] | 132 | { |
Gael Guennebaud | 6d1f5db | 2013-07-10 23:48:26 +0200 | [diff] [blame] | 133 | Index i = internal::random<Index>(0,rows-1); |
| 134 | Index j = internal::random<Index>(0,cols-1); |
Gael Guennebaud | 7706baf | 2011-09-08 13:42:54 +0200 | [diff] [blame] | 135 | if (m1.coeff(i,j)==Scalar(0)) |
| 136 | m2.insert(i,j) = m1(i,j) = internal::random<Scalar>(); |
| 137 | if(mode==3) |
| 138 | m2.reserve(r); |
| 139 | } |
Gael Guennebaud | 4ca89f3 | 2011-12-02 19:00:16 +0100 | [diff] [blame] | 140 | if(internal::random<int>()%2) |
| 141 | m2.makeCompressed(); |
Gael Guennebaud | 7706baf | 2011-09-08 13:42:54 +0200 | [diff] [blame] | 142 | VERIFY_IS_APPROX(m2,m1); |
| 143 | } |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 144 | |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 145 | // test basic computations |
| 146 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 147 | DenseMatrix refM1 = DenseMatrix::Zero(rows, cols); |
| 148 | DenseMatrix refM2 = DenseMatrix::Zero(rows, cols); |
| 149 | DenseMatrix refM3 = DenseMatrix::Zero(rows, cols); |
| 150 | DenseMatrix refM4 = DenseMatrix::Zero(rows, cols); |
| 151 | SparseMatrixType m1(rows, cols); |
| 152 | SparseMatrixType m2(rows, cols); |
| 153 | SparseMatrixType m3(rows, cols); |
| 154 | SparseMatrixType m4(rows, cols); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 155 | initSparse<Scalar>(density, refM1, m1); |
| 156 | initSparse<Scalar>(density, refM2, m2); |
| 157 | initSparse<Scalar>(density, refM3, m3); |
| 158 | initSparse<Scalar>(density, refM4, m4); |
| 159 | |
Gael Guennebaud | 4aac872 | 2014-07-22 12:54:03 +0200 | [diff] [blame] | 160 | VERIFY_IS_APPROX(m1*s1, refM1*s1); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 161 | VERIFY_IS_APPROX(m1+m2, refM1+refM2); |
| 162 | VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3); |
| 163 | VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2)); |
| 164 | VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2); |
| 165 | |
| 166 | VERIFY_IS_APPROX(m1*=s1, refM1*=s1); |
| 167 | VERIFY_IS_APPROX(m1/=s1, refM1/=s1); |
| 168 | |
| 169 | VERIFY_IS_APPROX(m1+=m2, refM1+=refM2); |
| 170 | VERIFY_IS_APPROX(m1-=m2, refM1-=refM2); |
| 171 | |
| 172 | if(SparseMatrixType::IsRowMajor) |
| 173 | VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0))); |
| 174 | else |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 175 | VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.col(0)), refM1.col(0).dot(refM2.col(0))); |
Gael Guennebaud | 3573a10 | 2014-02-17 13:46:17 +0100 | [diff] [blame] | 176 | |
| 177 | DenseVector rv = DenseVector::Random(m1.cols()); |
| 178 | DenseVector cv = DenseVector::Random(m1.rows()); |
| 179 | Index r = internal::random<Index>(0,m1.rows()-2); |
| 180 | Index c = internal::random<Index>(0,m1.cols()-1); |
| 181 | VERIFY_IS_APPROX(( m1.template block<1,Dynamic>(r,0,1,m1.cols()).dot(rv)) , refM1.row(r).dot(rv)); |
| 182 | VERIFY_IS_APPROX(m1.row(r).dot(rv), refM1.row(r).dot(rv)); |
| 183 | VERIFY_IS_APPROX(m1.col(c).dot(cv), refM1.col(c).dot(cv)); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 184 | |
| 185 | VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate()); |
| 186 | VERIFY_IS_APPROX(m1.real(), refM1.real()); |
| 187 | |
| 188 | refM4.setRandom(); |
| 189 | // sparse cwise* dense |
| 190 | VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4)); |
Gael Guennebaud | 9027508 | 2015-11-04 17:42:07 +0100 | [diff] [blame] | 191 | // dense cwise* sparse |
| 192 | VERIFY_IS_APPROX(refM4.cwiseProduct(m3), refM4.cwiseProduct(refM3)); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 193 | // VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4); |
| 194 | |
Gael Guennebaud | 15084cf | 2016-01-29 22:09:45 +0100 | [diff] [blame] | 195 | VERIFY_IS_APPROX(refM4 + m3, refM4 + refM3); |
| 196 | VERIFY_IS_APPROX(m3 + refM4, refM3 + refM4); |
| 197 | VERIFY_IS_APPROX(refM4 - m3, refM4 - refM3); |
| 198 | VERIFY_IS_APPROX(m3 - refM4, refM3 - refM4); |
| 199 | |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 200 | // test aliasing |
| 201 | VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1)); |
| 202 | VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval())); |
| 203 | VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval())); |
| 204 | VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1)); |
| 205 | } |
| 206 | |
| 207 | // test transpose |
| 208 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 209 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols); |
| 210 | SparseMatrixType m2(rows, cols); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 211 | initSparse<Scalar>(density, refMat2, m2); |
| 212 | VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval()); |
| 213 | VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose()); |
| 214 | |
| 215 | VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint()); |
Gael Guennebaud | ff46ec0 | 2014-09-22 23:33:28 +0200 | [diff] [blame] | 216 | |
| 217 | // check isApprox handles opposite storage order |
| 218 | typename Transpose<SparseMatrixType>::PlainObject m3(m2); |
| 219 | VERIFY(m2.isApprox(m3)); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 220 | } |
| 221 | |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 222 | // test prune |
| 223 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 224 | SparseMatrixType m2(rows, cols); |
| 225 | DenseMatrix refM2(rows, cols); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 226 | refM2.setZero(); |
| 227 | int countFalseNonZero = 0; |
| 228 | int countTrueNonZero = 0; |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 229 | m2.reserve(VectorXi::Constant(m2.outerSize(), int(m2.innerSize()))); |
| 230 | for (Index j=0; j<m2.cols(); ++j) |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 231 | { |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 232 | for (Index i=0; i<m2.rows(); ++i) |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 233 | { |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 234 | float x = internal::random<float>(0,1); |
Christoph Hertzberg | dacb469 | 2016-05-05 13:35:45 +0200 | [diff] [blame] | 235 | if (x<0.1f) |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 236 | { |
| 237 | // do nothing |
| 238 | } |
Christoph Hertzberg | dacb469 | 2016-05-05 13:35:45 +0200 | [diff] [blame] | 239 | else if (x<0.5f) |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 240 | { |
| 241 | countFalseNonZero++; |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 242 | m2.insert(i,j) = Scalar(0); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 243 | } |
| 244 | else |
| 245 | { |
| 246 | countTrueNonZero++; |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 247 | m2.insert(i,j) = Scalar(1); |
| 248 | refM2(i,j) = Scalar(1); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 249 | } |
| 250 | } |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 251 | } |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 252 | if(internal::random<bool>()) |
| 253 | m2.makeCompressed(); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 254 | VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros()); |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 255 | if(countTrueNonZero>0) |
| 256 | VERIFY_IS_APPROX(m2, refM2); |
Hauke Heibel | d204ec4 | 2010-11-02 14:33:33 +0100 | [diff] [blame] | 257 | m2.prune(Scalar(1)); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 258 | VERIFY(countTrueNonZero==m2.nonZeros()); |
| 259 | VERIFY_IS_APPROX(m2, refM2); |
| 260 | } |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 261 | |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 262 | // test setFromTriplets |
| 263 | { |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 264 | typedef Triplet<Scalar,StorageIndex> TripletType; |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 265 | std::vector<TripletType> triplets; |
Gael Guennebaud | b47ef14 | 2014-07-08 16:47:11 +0200 | [diff] [blame] | 266 | Index ntriplets = rows*cols; |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 267 | triplets.reserve(ntriplets); |
Gael Guennebaud | b4c79ee | 2015-10-13 11:30:41 +0200 | [diff] [blame] | 268 | DenseMatrix refMat_sum = DenseMatrix::Zero(rows,cols); |
| 269 | DenseMatrix refMat_prod = DenseMatrix::Zero(rows,cols); |
| 270 | DenseMatrix refMat_last = DenseMatrix::Zero(rows,cols); |
| 271 | |
Gael Guennebaud | b47ef14 | 2014-07-08 16:47:11 +0200 | [diff] [blame] | 272 | for(Index i=0;i<ntriplets;++i) |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 273 | { |
Gael Guennebaud | aa6c516 | 2015-02-16 13:19:05 +0100 | [diff] [blame] | 274 | StorageIndex r = internal::random<StorageIndex>(0,StorageIndex(rows-1)); |
| 275 | StorageIndex c = internal::random<StorageIndex>(0,StorageIndex(cols-1)); |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 276 | Scalar v = internal::random<Scalar>(); |
Gael Guennebaud | 18e3ac0 | 2012-01-31 09:14:01 +0100 | [diff] [blame] | 277 | triplets.push_back(TripletType(r,c,v)); |
Gael Guennebaud | b4c79ee | 2015-10-13 11:30:41 +0200 | [diff] [blame] | 278 | refMat_sum(r,c) += v; |
| 279 | if(std::abs(refMat_prod(r,c))==0) |
| 280 | refMat_prod(r,c) = v; |
| 281 | else |
| 282 | refMat_prod(r,c) *= v; |
| 283 | refMat_last(r,c) = v; |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 284 | } |
| 285 | SparseMatrixType m(rows,cols); |
| 286 | m.setFromTriplets(triplets.begin(), triplets.end()); |
Gael Guennebaud | b4c79ee | 2015-10-13 11:30:41 +0200 | [diff] [blame] | 287 | VERIFY_IS_APPROX(m, refMat_sum); |
| 288 | |
| 289 | m.setFromTriplets(triplets.begin(), triplets.end(), std::multiplies<Scalar>()); |
| 290 | VERIFY_IS_APPROX(m, refMat_prod); |
| 291 | #if (defined(__cplusplus) && __cplusplus >= 201103L) |
| 292 | m.setFromTriplets(triplets.begin(), triplets.end(), [] (Scalar,Scalar b) { return b; }); |
| 293 | VERIFY_IS_APPROX(m, refMat_last); |
| 294 | #endif |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 295 | } |
Gael Guennebaud | 3af29ca | 2015-02-09 10:23:45 +0100 | [diff] [blame] | 296 | |
| 297 | // test Map |
| 298 | { |
| 299 | DenseMatrix refMat2(rows, cols), refMat3(rows, cols); |
| 300 | SparseMatrixType m2(rows, cols), m3(rows, cols); |
| 301 | initSparse<Scalar>(density, refMat2, m2); |
| 302 | initSparse<Scalar>(density, refMat3, m3); |
| 303 | { |
| 304 | Map<SparseMatrixType> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr()); |
| 305 | Map<SparseMatrixType> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr()); |
| 306 | VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3); |
| 307 | VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3); |
| 308 | } |
| 309 | { |
Gael Guennebaud | fe51319 | 2015-02-13 10:03:53 +0100 | [diff] [blame] | 310 | MappedSparseMatrix<Scalar,SparseMatrixType::Options,StorageIndex> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr()); |
| 311 | MappedSparseMatrix<Scalar,SparseMatrixType::Options,StorageIndex> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr()); |
Gael Guennebaud | 3af29ca | 2015-02-09 10:23:45 +0100 | [diff] [blame] | 312 | VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3); |
| 313 | VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3); |
| 314 | } |
| 315 | } |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 316 | |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 317 | // test triangularView |
| 318 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 319 | DenseMatrix refMat2(rows, cols), refMat3(rows, cols); |
| 320 | SparseMatrixType m2(rows, cols), m3(rows, cols); |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 321 | initSparse<Scalar>(density, refMat2, m2); |
| 322 | refMat3 = refMat2.template triangularView<Lower>(); |
| 323 | m3 = m2.template triangularView<Lower>(); |
| 324 | VERIFY_IS_APPROX(m3, refMat3); |
| 325 | |
| 326 | refMat3 = refMat2.template triangularView<Upper>(); |
| 327 | m3 = m2.template triangularView<Upper>(); |
| 328 | VERIFY_IS_APPROX(m3, refMat3); |
| 329 | |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 330 | { |
| 331 | refMat3 = refMat2.template triangularView<UnitUpper>(); |
| 332 | m3 = m2.template triangularView<UnitUpper>(); |
| 333 | VERIFY_IS_APPROX(m3, refMat3); |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 334 | |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 335 | refMat3 = refMat2.template triangularView<UnitLower>(); |
| 336 | m3 = m2.template triangularView<UnitLower>(); |
| 337 | VERIFY_IS_APPROX(m3, refMat3); |
| 338 | } |
Desire NUENTSA | 4cd8245 | 2013-06-11 14:42:29 +0200 | [diff] [blame] | 339 | |
| 340 | refMat3 = refMat2.template triangularView<StrictlyUpper>(); |
| 341 | m3 = m2.template triangularView<StrictlyUpper>(); |
| 342 | VERIFY_IS_APPROX(m3, refMat3); |
| 343 | |
| 344 | refMat3 = refMat2.template triangularView<StrictlyLower>(); |
| 345 | m3 = m2.template triangularView<StrictlyLower>(); |
| 346 | VERIFY_IS_APPROX(m3, refMat3); |
Gael Guennebaud | f6b1dee | 2015-11-04 17:02:32 +0100 | [diff] [blame] | 347 | |
| 348 | // check sparse-traingular to dense |
| 349 | refMat3 = m2.template triangularView<StrictlyUpper>(); |
| 350 | VERIFY_IS_APPROX(refMat3, DenseMatrix(refMat2.template triangularView<StrictlyUpper>())); |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 351 | } |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 352 | |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 353 | // test selfadjointView |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 354 | if(!SparseMatrixType::IsRowMajor) |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 355 | { |
| 356 | DenseMatrix refMat2(rows, rows), refMat3(rows, rows); |
| 357 | SparseMatrixType m2(rows, rows), m3(rows, rows); |
| 358 | initSparse<Scalar>(density, refMat2, m2); |
| 359 | refMat3 = refMat2.template selfadjointView<Lower>(); |
| 360 | m3 = m2.template selfadjointView<Lower>(); |
| 361 | VERIFY_IS_APPROX(m3, refMat3); |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 362 | |
| 363 | // selfadjointView only works for square matrices: |
| 364 | SparseMatrixType m4(rows, rows+1); |
| 365 | VERIFY_RAISES_ASSERT(m4.template selfadjointView<Lower>()); |
| 366 | VERIFY_RAISES_ASSERT(m4.template selfadjointView<Upper>()); |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 367 | } |
| 368 | |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 369 | // test sparseView |
| 370 | { |
| 371 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows); |
| 372 | SparseMatrixType m2(rows, rows); |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 373 | initSparse<Scalar>(density, refMat2, m2); |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 374 | VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval()); |
Gael Guennebaud | 8456bbb | 2016-05-18 16:53:28 +0200 | [diff] [blame^] | 375 | |
| 376 | // sparse view on expressions: |
| 377 | VERIFY_IS_APPROX((s1*m2).eval(), (s1*refMat2).sparseView().eval()); |
| 378 | VERIFY_IS_APPROX((m2+m2).eval(), (refMat2+refMat2).sparseView().eval()); |
| 379 | VERIFY_IS_APPROX((m2*m2).eval(), (refMat2.lazyProduct(refMat2)).sparseView().eval()); |
| 380 | VERIFY_IS_APPROX((m2*m2).eval(), (refMat2*refMat2).sparseView().eval()); |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 381 | } |
Gael Guennebaud | 82f9aa1 | 2011-12-04 21:49:21 +0100 | [diff] [blame] | 382 | |
| 383 | // test diagonal |
| 384 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 385 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols); |
| 386 | SparseMatrixType m2(rows, cols); |
Gael Guennebaud | 82f9aa1 | 2011-12-04 21:49:21 +0100 | [diff] [blame] | 387 | initSparse<Scalar>(density, refMat2, m2); |
| 388 | VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval()); |
Gael Guennebaud | b26e697 | 2014-12-01 14:41:39 +0100 | [diff] [blame] | 389 | VERIFY_IS_APPROX(const_cast<const SparseMatrixType&>(m2).diagonal(), refMat2.diagonal().eval()); |
| 390 | |
| 391 | initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag); |
| 392 | m2.diagonal() += refMat2.diagonal(); |
| 393 | refMat2.diagonal() += refMat2.diagonal(); |
| 394 | VERIFY_IS_APPROX(m2, refMat2); |
Gael Guennebaud | 82f9aa1 | 2011-12-04 21:49:21 +0100 | [diff] [blame] | 395 | } |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 396 | |
Gael Guennebaud | 62f21e2 | 2015-06-24 17:55:00 +0200 | [diff] [blame] | 397 | // test diagonal to sparse |
| 398 | { |
| 399 | DenseVector d = DenseVector::Random(rows); |
| 400 | DenseMatrix refMat2 = d.asDiagonal(); |
| 401 | SparseMatrixType m2(rows, rows); |
| 402 | m2 = d.asDiagonal(); |
| 403 | VERIFY_IS_APPROX(m2, refMat2); |
Gael Guennebaud | 4c8cd13 | 2015-06-24 18:11:06 +0200 | [diff] [blame] | 404 | SparseMatrixType m3(d.asDiagonal()); |
| 405 | VERIFY_IS_APPROX(m3, refMat2); |
Gael Guennebaud | 62f21e2 | 2015-06-24 17:55:00 +0200 | [diff] [blame] | 406 | refMat2 += d.asDiagonal(); |
| 407 | m2 += d.asDiagonal(); |
| 408 | VERIFY_IS_APPROX(m2, refMat2); |
| 409 | } |
| 410 | |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 411 | // test conservative resize |
| 412 | { |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 413 | std::vector< std::pair<StorageIndex,StorageIndex> > inc; |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 414 | if(rows > 3 && cols > 2) |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 415 | inc.push_back(std::pair<StorageIndex,StorageIndex>(-3,-2)); |
| 416 | inc.push_back(std::pair<StorageIndex,StorageIndex>(0,0)); |
| 417 | inc.push_back(std::pair<StorageIndex,StorageIndex>(3,2)); |
| 418 | inc.push_back(std::pair<StorageIndex,StorageIndex>(3,0)); |
| 419 | inc.push_back(std::pair<StorageIndex,StorageIndex>(0,3)); |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 420 | |
| 421 | for(size_t i = 0; i< inc.size(); i++) { |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 422 | StorageIndex incRows = inc[i].first; |
| 423 | StorageIndex incCols = inc[i].second; |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 424 | SparseMatrixType m1(rows, cols); |
| 425 | DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols); |
| 426 | initSparse<Scalar>(density, refMat1, m1); |
| 427 | |
| 428 | m1.conservativeResize(rows+incRows, cols+incCols); |
| 429 | refMat1.conservativeResize(rows+incRows, cols+incCols); |
| 430 | if (incRows > 0) refMat1.bottomRows(incRows).setZero(); |
| 431 | if (incCols > 0) refMat1.rightCols(incCols).setZero(); |
| 432 | |
| 433 | VERIFY_IS_APPROX(m1, refMat1); |
| 434 | |
| 435 | // Insert new values |
| 436 | if (incRows > 0) |
Gael Guennebaud | 7ee378d | 2013-07-12 16:40:02 +0200 | [diff] [blame] | 437 | m1.insert(m1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1; |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 438 | if (incCols > 0) |
Gael Guennebaud | 7ee378d | 2013-07-12 16:40:02 +0200 | [diff] [blame] | 439 | m1.insert(0, m1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1; |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 440 | |
| 441 | VERIFY_IS_APPROX(m1, refMat1); |
| 442 | |
| 443 | |
| 444 | } |
| 445 | } |
Desire NUENTSA | 4cd8245 | 2013-06-11 14:42:29 +0200 | [diff] [blame] | 446 | |
| 447 | // test Identity matrix |
| 448 | { |
| 449 | DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows); |
| 450 | SparseMatrixType m1(rows, rows); |
| 451 | m1.setIdentity(); |
| 452 | VERIFY_IS_APPROX(m1, refMat1); |
Gael Guennebaud | 8a211bb | 2015-10-25 22:01:58 +0100 | [diff] [blame] | 453 | for(int k=0; k<rows*rows/4; ++k) |
| 454 | { |
| 455 | Index i = internal::random<Index>(0,rows-1); |
| 456 | Index j = internal::random<Index>(0,rows-1); |
Gael Guennebaud | 73f692d | 2015-10-27 11:01:37 +0100 | [diff] [blame] | 457 | Scalar v = internal::random<Scalar>(); |
Gael Guennebaud | 8a211bb | 2015-10-25 22:01:58 +0100 | [diff] [blame] | 458 | m1.coeffRef(i,j) = v; |
| 459 | refMat1.coeffRef(i,j) = v; |
| 460 | VERIFY_IS_APPROX(m1, refMat1); |
| 461 | if(internal::random<Index>(0,10)<2) |
| 462 | m1.makeCompressed(); |
| 463 | } |
| 464 | m1.setIdentity(); |
| 465 | refMat1.setIdentity(); |
| 466 | VERIFY_IS_APPROX(m1, refMat1); |
Desire NUENTSA | 4cd8245 | 2013-06-11 14:42:29 +0200 | [diff] [blame] | 467 | } |
Gael Guennebaud | ec46970 | 2016-02-01 15:04:33 +0100 | [diff] [blame] | 468 | |
| 469 | // test array/vector of InnerIterator |
| 470 | { |
| 471 | typedef typename SparseMatrixType::InnerIterator IteratorType; |
| 472 | |
| 473 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols); |
| 474 | SparseMatrixType m2(rows, cols); |
| 475 | initSparse<Scalar>(density, refMat2, m2); |
| 476 | IteratorType static_array[2]; |
| 477 | static_array[0] = IteratorType(m2,0); |
| 478 | static_array[1] = IteratorType(m2,m2.outerSize()-1); |
| 479 | VERIFY( static_array[0] || m2.innerVector(static_array[0].outer()).nonZeros() == 0 ); |
| 480 | VERIFY( static_array[1] || m2.innerVector(static_array[1].outer()).nonZeros() == 0 ); |
| 481 | if(static_array[0] && static_array[1]) |
| 482 | { |
| 483 | ++(static_array[1]); |
| 484 | static_array[1] = IteratorType(m2,0); |
| 485 | VERIFY( static_array[1] ); |
| 486 | VERIFY( static_array[1].index() == static_array[0].index() ); |
| 487 | VERIFY( static_array[1].outer() == static_array[0].outer() ); |
| 488 | VERIFY( static_array[1].value() == static_array[0].value() ); |
| 489 | } |
| 490 | |
| 491 | std::vector<IteratorType> iters(2); |
| 492 | iters[0] = IteratorType(m2,0); |
| 493 | iters[1] = IteratorType(m2,m2.outerSize()-1); |
| 494 | } |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 495 | } |
| 496 | |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 497 | |
Christoph Hertzberg | c5a3777 | 2014-10-31 17:19:05 +0100 | [diff] [blame] | 498 | template<typename SparseMatrixType> |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 499 | void big_sparse_triplet(Index rows, Index cols, double density) { |
| 500 | typedef typename SparseMatrixType::StorageIndex StorageIndex; |
| 501 | typedef typename SparseMatrixType::Scalar Scalar; |
| 502 | typedef Triplet<Scalar,Index> TripletType; |
| 503 | std::vector<TripletType> triplets; |
| 504 | double nelements = density * rows*cols; |
| 505 | VERIFY(nelements>=0 && nelements < NumTraits<StorageIndex>::highest()); |
| 506 | Index ntriplets = Index(nelements); |
| 507 | triplets.reserve(ntriplets); |
| 508 | Scalar sum = Scalar(0); |
| 509 | for(Index i=0;i<ntriplets;++i) |
| 510 | { |
| 511 | Index r = internal::random<Index>(0,rows-1); |
| 512 | Index c = internal::random<Index>(0,cols-1); |
| 513 | Scalar v = internal::random<Scalar>(); |
| 514 | triplets.push_back(TripletType(r,c,v)); |
| 515 | sum += v; |
| 516 | } |
| 517 | SparseMatrixType m(rows,cols); |
| 518 | m.setFromTriplets(triplets.begin(), triplets.end()); |
| 519 | VERIFY(m.nonZeros() <= ntriplets); |
| 520 | VERIFY_IS_APPROX(sum, m.sum()); |
Christoph Hertzberg | c5a3777 | 2014-10-31 17:19:05 +0100 | [diff] [blame] | 521 | } |
| 522 | |
| 523 | |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 524 | void test_sparse_basic() |
| 525 | { |
| 526 | for(int i = 0; i < g_repeat; i++) { |
Gael Guennebaud | c43154b | 2015-03-04 10:16:46 +0100 | [diff] [blame] | 527 | int r = Eigen::internal::random<int>(1,200), c = Eigen::internal::random<int>(1,200); |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 528 | if(Eigen::internal::random<int>(0,4) == 0) { |
| 529 | r = c; // check square matrices in 25% of tries |
| 530 | } |
| 531 | EIGEN_UNUSED_VARIABLE(r+c); |
| 532 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(1, 1)) )); |
Gael Guennebaud | 91fe150 | 2011-06-07 11:28:16 +0200 | [diff] [blame] | 533 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) )); |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 534 | CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(r, c)) )); |
| 535 | CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(r, c)) )); |
| 536 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(r, c)) )); |
Gael Guennebaud | d7698c1 | 2015-03-19 15:11:05 +0100 | [diff] [blame] | 537 | CALL_SUBTEST_5(( sparse_basic(SparseMatrix<double,ColMajor,long int>(r, c)) )); |
| 538 | CALL_SUBTEST_5(( sparse_basic(SparseMatrix<double,RowMajor,long int>(r, c)) )); |
Gael Guennebaud | 47e8902 | 2013-06-10 10:34:03 +0200 | [diff] [blame] | 539 | |
Gael Guennebaud | c43154b | 2015-03-04 10:16:46 +0100 | [diff] [blame] | 540 | r = Eigen::internal::random<int>(1,100); |
| 541 | c = Eigen::internal::random<int>(1,100); |
| 542 | if(Eigen::internal::random<int>(0,4) == 0) { |
| 543 | r = c; // check square matrices in 25% of tries |
| 544 | } |
| 545 | |
Gael Guennebaud | d7698c1 | 2015-03-19 15:11:05 +0100 | [diff] [blame] | 546 | CALL_SUBTEST_6(( sparse_basic(SparseMatrix<double,ColMajor,short int>(short(r), short(c))) )); |
| 547 | CALL_SUBTEST_6(( sparse_basic(SparseMatrix<double,RowMajor,short int>(short(r), short(c))) )); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 548 | } |
Christoph Hertzberg | c5a3777 | 2014-10-31 17:19:05 +0100 | [diff] [blame] | 549 | |
| 550 | // Regression test for bug 900: (manually insert higher values here, if you have enough RAM): |
| 551 | CALL_SUBTEST_3((big_sparse_triplet<SparseMatrix<float, RowMajor, int> >(10000, 10000, 0.125))); |
| 552 | CALL_SUBTEST_4((big_sparse_triplet<SparseMatrix<double, ColMajor, long int> >(10000, 10000, 0.125))); |
Gael Guennebaud | bfd6ee6 | 2015-11-06 15:05:37 +0100 | [diff] [blame] | 553 | |
| 554 | // Regression test for bug 1105 |
Christoph Hertzberg | 7268b10 | 2016-05-11 19:41:53 +0200 | [diff] [blame] | 555 | #ifdef EIGEN_TEST_PART_7 |
Gael Guennebaud | bfd6ee6 | 2015-11-06 15:05:37 +0100 | [diff] [blame] | 556 | { |
| 557 | int n = Eigen::internal::random<int>(200,600); |
| 558 | SparseMatrix<std::complex<double>,0, long> mat(n, n); |
| 559 | std::complex<double> val; |
| 560 | |
| 561 | for(int i=0; i<n; ++i) |
| 562 | { |
| 563 | mat.coeffRef(i, i%(n/10)) = val; |
| 564 | VERIFY(mat.data().allocatedSize()<20*n); |
| 565 | } |
| 566 | } |
| 567 | #endif |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 568 | } |