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 | 2c1b56f | 2016-05-31 10:56:53 +0200 | [diff] [blame] | 160 | if(internal::random<bool>()) |
| 161 | m1.makeCompressed(); |
| 162 | |
Gael Guennebaud | 4aac872 | 2014-07-22 12:54:03 +0200 | [diff] [blame] | 163 | VERIFY_IS_APPROX(m1*s1, refM1*s1); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 164 | VERIFY_IS_APPROX(m1+m2, refM1+refM2); |
| 165 | VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3); |
| 166 | VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2)); |
| 167 | VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2); |
| 168 | |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 169 | if(SparseMatrixType::IsRowMajor) |
| 170 | VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0))); |
| 171 | else |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 172 | 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] | 173 | |
| 174 | DenseVector rv = DenseVector::Random(m1.cols()); |
| 175 | DenseVector cv = DenseVector::Random(m1.rows()); |
| 176 | Index r = internal::random<Index>(0,m1.rows()-2); |
| 177 | Index c = internal::random<Index>(0,m1.cols()-1); |
| 178 | VERIFY_IS_APPROX(( m1.template block<1,Dynamic>(r,0,1,m1.cols()).dot(rv)) , refM1.row(r).dot(rv)); |
| 179 | VERIFY_IS_APPROX(m1.row(r).dot(rv), refM1.row(r).dot(rv)); |
| 180 | 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] | 181 | |
| 182 | VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate()); |
| 183 | VERIFY_IS_APPROX(m1.real(), refM1.real()); |
| 184 | |
| 185 | refM4.setRandom(); |
| 186 | // sparse cwise* dense |
| 187 | VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4)); |
Gael Guennebaud | 9027508 | 2015-11-04 17:42:07 +0100 | [diff] [blame] | 188 | // dense cwise* sparse |
| 189 | VERIFY_IS_APPROX(refM4.cwiseProduct(m3), refM4.cwiseProduct(refM3)); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 190 | // VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4); |
| 191 | |
Gael Guennebaud | 15084cf | 2016-01-29 22:09:45 +0100 | [diff] [blame] | 192 | VERIFY_IS_APPROX(refM4 + m3, refM4 + refM3); |
| 193 | VERIFY_IS_APPROX(m3 + refM4, refM3 + refM4); |
| 194 | VERIFY_IS_APPROX(refM4 - m3, refM4 - refM3); |
| 195 | VERIFY_IS_APPROX(m3 - refM4, refM3 - refM4); |
| 196 | |
Gael Guennebaud | 2c1b56f | 2016-05-31 10:56:53 +0200 | [diff] [blame] | 197 | VERIFY_IS_APPROX(m1.sum(), refM1.sum()); |
| 198 | |
| 199 | VERIFY_IS_APPROX(m1*=s1, refM1*=s1); |
| 200 | VERIFY_IS_APPROX(m1/=s1, refM1/=s1); |
| 201 | |
| 202 | VERIFY_IS_APPROX(m1+=m2, refM1+=refM2); |
| 203 | VERIFY_IS_APPROX(m1-=m2, refM1-=refM2); |
| 204 | |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 205 | // test aliasing |
| 206 | VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1)); |
| 207 | VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval())); |
| 208 | VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval())); |
| 209 | VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1)); |
Gael Guennebaud | 7e029d1 | 2016-08-29 12:06:37 +0200 | [diff] [blame] | 210 | |
| 211 | if(m1.isCompressed()) |
| 212 | { |
| 213 | VERIFY_IS_APPROX(m1.coeffs().sum(), m1.sum()); |
| 214 | m1.coeffs() += s1; |
| 215 | for(Index j = 0; j<m1.outerSize(); ++j) |
| 216 | for(typename SparseMatrixType::InnerIterator it(m1,j); it; ++it) |
| 217 | refM1(it.row(), it.col()) += s1; |
| 218 | VERIFY_IS_APPROX(m1, refM1); |
| 219 | } |
Gael Guennebaud | 2e334f5 | 2016-11-14 18:47:02 +0100 | [diff] [blame] | 220 | |
| 221 | // and/or |
| 222 | { |
| 223 | typedef SparseMatrix<bool, SparseMatrixType::Options, typename SparseMatrixType::StorageIndex> SpBool; |
| 224 | SpBool mb1 = m1.real().template cast<bool>(); |
| 225 | SpBool mb2 = m2.real().template cast<bool>(); |
| 226 | VERIFY_IS_EQUAL(mb1.template cast<int>().sum(), refM1.real().template cast<bool>().count()); |
| 227 | VERIFY_IS_EQUAL((mb1 && mb2).template cast<int>().sum(), (refM1.real().template cast<bool>() && refM2.real().template cast<bool>()).count()); |
| 228 | VERIFY_IS_EQUAL((mb1 || mb2).template cast<int>().sum(), (refM1.real().template cast<bool>() || refM2.real().template cast<bool>()).count()); |
| 229 | SpBool mb3 = mb1 && mb2; |
| 230 | if(mb1.coeffs().all() && mb2.coeffs().all()) |
| 231 | { |
| 232 | VERIFY_IS_EQUAL(mb3.nonZeros(), (refM1.real().template cast<bool>() && refM2.real().template cast<bool>()).count()); |
| 233 | } |
| 234 | } |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 235 | } |
| 236 | |
Gael Guennebaud | eedb87f | 2016-11-14 14:05:53 +0100 | [diff] [blame] | 237 | // test reverse iterators |
| 238 | { |
| 239 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols); |
| 240 | SparseMatrixType m2(rows, cols); |
| 241 | initSparse<Scalar>(density, refMat2, m2); |
| 242 | std::vector<Scalar> ref_value(m2.innerSize()); |
| 243 | std::vector<Index> ref_index(m2.innerSize()); |
| 244 | if(internal::random<bool>()) |
| 245 | m2.makeCompressed(); |
| 246 | for(Index j = 0; j<m2.outerSize(); ++j) |
| 247 | { |
| 248 | Index count_forward = 0; |
| 249 | |
| 250 | for(typename SparseMatrixType::InnerIterator it(m2,j); it; ++it) |
| 251 | { |
| 252 | ref_value[ref_value.size()-1-count_forward] = it.value(); |
| 253 | ref_index[ref_index.size()-1-count_forward] = it.index(); |
| 254 | count_forward++; |
| 255 | } |
| 256 | Index count_reverse = 0; |
| 257 | for(typename SparseMatrixType::ReverseInnerIterator it(m2,j); it; --it) |
| 258 | { |
| 259 | VERIFY_IS_APPROX( std::abs(ref_value[ref_value.size()-count_forward+count_reverse])+1, std::abs(it.value())+1); |
| 260 | VERIFY_IS_EQUAL( ref_index[ref_index.size()-count_forward+count_reverse] , it.index()); |
| 261 | count_reverse++; |
| 262 | } |
| 263 | VERIFY_IS_EQUAL(count_forward, count_reverse); |
| 264 | } |
| 265 | } |
| 266 | |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 267 | // test transpose |
| 268 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 269 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols); |
| 270 | SparseMatrixType m2(rows, cols); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 271 | initSparse<Scalar>(density, refMat2, m2); |
| 272 | VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval()); |
| 273 | VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose()); |
| 274 | |
| 275 | VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint()); |
Gael Guennebaud | ff46ec0 | 2014-09-22 23:33:28 +0200 | [diff] [blame] | 276 | |
| 277 | // check isApprox handles opposite storage order |
| 278 | typename Transpose<SparseMatrixType>::PlainObject m3(m2); |
| 279 | VERIFY(m2.isApprox(m3)); |
Gael Guennebaud | 4e60283 | 2012-11-16 09:02:50 +0100 | [diff] [blame] | 280 | } |
| 281 | |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 282 | // test prune |
| 283 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 284 | SparseMatrixType m2(rows, cols); |
| 285 | DenseMatrix refM2(rows, cols); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 286 | refM2.setZero(); |
| 287 | int countFalseNonZero = 0; |
| 288 | int countTrueNonZero = 0; |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 289 | m2.reserve(VectorXi::Constant(m2.outerSize(), int(m2.innerSize()))); |
| 290 | for (Index j=0; j<m2.cols(); ++j) |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 291 | { |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 292 | for (Index i=0; i<m2.rows(); ++i) |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 293 | { |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 294 | float x = internal::random<float>(0,1); |
Christoph Hertzberg | dacb469 | 2016-05-05 13:35:45 +0200 | [diff] [blame] | 295 | if (x<0.1f) |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 296 | { |
| 297 | // do nothing |
| 298 | } |
Christoph Hertzberg | dacb469 | 2016-05-05 13:35:45 +0200 | [diff] [blame] | 299 | else if (x<0.5f) |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 300 | { |
| 301 | countFalseNonZero++; |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 302 | m2.insert(i,j) = Scalar(0); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 303 | } |
| 304 | else |
| 305 | { |
| 306 | countTrueNonZero++; |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 307 | m2.insert(i,j) = Scalar(1); |
| 308 | refM2(i,j) = Scalar(1); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 309 | } |
| 310 | } |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 311 | } |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 312 | if(internal::random<bool>()) |
| 313 | m2.makeCompressed(); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 314 | VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros()); |
Gael Guennebaud | a44d91a | 2015-10-13 10:53:38 +0200 | [diff] [blame] | 315 | if(countTrueNonZero>0) |
| 316 | VERIFY_IS_APPROX(m2, refM2); |
Hauke Heibel | d204ec4 | 2010-11-02 14:33:33 +0100 | [diff] [blame] | 317 | m2.prune(Scalar(1)); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 318 | VERIFY(countTrueNonZero==m2.nonZeros()); |
| 319 | VERIFY_IS_APPROX(m2, refM2); |
| 320 | } |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 321 | |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 322 | // test setFromTriplets |
| 323 | { |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 324 | typedef Triplet<Scalar,StorageIndex> TripletType; |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 325 | std::vector<TripletType> triplets; |
Gael Guennebaud | b47ef14 | 2014-07-08 16:47:11 +0200 | [diff] [blame] | 326 | Index ntriplets = rows*cols; |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 327 | triplets.reserve(ntriplets); |
Gael Guennebaud | b4c79ee | 2015-10-13 11:30:41 +0200 | [diff] [blame] | 328 | DenseMatrix refMat_sum = DenseMatrix::Zero(rows,cols); |
| 329 | DenseMatrix refMat_prod = DenseMatrix::Zero(rows,cols); |
| 330 | DenseMatrix refMat_last = DenseMatrix::Zero(rows,cols); |
| 331 | |
Gael Guennebaud | b47ef14 | 2014-07-08 16:47:11 +0200 | [diff] [blame] | 332 | for(Index i=0;i<ntriplets;++i) |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 333 | { |
Gael Guennebaud | aa6c516 | 2015-02-16 13:19:05 +0100 | [diff] [blame] | 334 | StorageIndex r = internal::random<StorageIndex>(0,StorageIndex(rows-1)); |
| 335 | StorageIndex c = internal::random<StorageIndex>(0,StorageIndex(cols-1)); |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 336 | Scalar v = internal::random<Scalar>(); |
Gael Guennebaud | 18e3ac0 | 2012-01-31 09:14:01 +0100 | [diff] [blame] | 337 | triplets.push_back(TripletType(r,c,v)); |
Gael Guennebaud | b4c79ee | 2015-10-13 11:30:41 +0200 | [diff] [blame] | 338 | refMat_sum(r,c) += v; |
| 339 | if(std::abs(refMat_prod(r,c))==0) |
| 340 | refMat_prod(r,c) = v; |
| 341 | else |
| 342 | refMat_prod(r,c) *= v; |
| 343 | refMat_last(r,c) = v; |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 344 | } |
| 345 | SparseMatrixType m(rows,cols); |
| 346 | m.setFromTriplets(triplets.begin(), triplets.end()); |
Gael Guennebaud | b4c79ee | 2015-10-13 11:30:41 +0200 | [diff] [blame] | 347 | VERIFY_IS_APPROX(m, refMat_sum); |
| 348 | |
| 349 | m.setFromTriplets(triplets.begin(), triplets.end(), std::multiplies<Scalar>()); |
| 350 | VERIFY_IS_APPROX(m, refMat_prod); |
| 351 | #if (defined(__cplusplus) && __cplusplus >= 201103L) |
| 352 | m.setFromTriplets(triplets.begin(), triplets.end(), [] (Scalar,Scalar b) { return b; }); |
| 353 | VERIFY_IS_APPROX(m, refMat_last); |
| 354 | #endif |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 355 | } |
Gael Guennebaud | 3af29ca | 2015-02-09 10:23:45 +0100 | [diff] [blame] | 356 | |
| 357 | // test Map |
| 358 | { |
| 359 | DenseMatrix refMat2(rows, cols), refMat3(rows, cols); |
| 360 | SparseMatrixType m2(rows, cols), m3(rows, cols); |
| 361 | initSparse<Scalar>(density, refMat2, m2); |
| 362 | initSparse<Scalar>(density, refMat3, m3); |
| 363 | { |
| 364 | Map<SparseMatrixType> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr()); |
| 365 | Map<SparseMatrixType> mapMat3(m3.rows(), m3.cols(), m3.nonZeros(), m3.outerIndexPtr(), m3.innerIndexPtr(), m3.valuePtr(), m3.innerNonZeroPtr()); |
| 366 | VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3); |
| 367 | VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3); |
| 368 | } |
| 369 | { |
Gael Guennebaud | fe51319 | 2015-02-13 10:03:53 +0100 | [diff] [blame] | 370 | MappedSparseMatrix<Scalar,SparseMatrixType::Options,StorageIndex> mapMat2(m2.rows(), m2.cols(), m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr()); |
| 371 | 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] | 372 | VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3); |
| 373 | VERIFY_IS_APPROX(mapMat2+mapMat3, refMat2+refMat3); |
| 374 | } |
Gael Guennebaud | 757971e | 2016-07-26 09:40:19 +0200 | [diff] [blame] | 375 | |
| 376 | Index i = internal::random<Index>(0,rows-1); |
| 377 | Index j = internal::random<Index>(0,cols-1); |
| 378 | m2.coeffRef(i,j) = 123; |
| 379 | if(internal::random<bool>()) |
| 380 | m2.makeCompressed(); |
| 381 | Map<SparseMatrixType> mapMat2(rows, cols, m2.nonZeros(), m2.outerIndexPtr(), m2.innerIndexPtr(), m2.valuePtr(), m2.innerNonZeroPtr()); |
| 382 | VERIFY_IS_EQUAL(m2.coeff(i,j),Scalar(123)); |
| 383 | VERIFY_IS_EQUAL(mapMat2.coeff(i,j),Scalar(123)); |
| 384 | mapMat2.coeffRef(i,j) = -123; |
| 385 | VERIFY_IS_EQUAL(m2.coeff(i,j),Scalar(-123)); |
Gael Guennebaud | 3af29ca | 2015-02-09 10:23:45 +0100 | [diff] [blame] | 386 | } |
Gael Guennebaud | 8713807 | 2012-01-28 11:13:59 +0100 | [diff] [blame] | 387 | |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 388 | // test triangularView |
| 389 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 390 | DenseMatrix refMat2(rows, cols), refMat3(rows, cols); |
| 391 | SparseMatrixType m2(rows, cols), m3(rows, cols); |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 392 | initSparse<Scalar>(density, refMat2, m2); |
| 393 | refMat3 = refMat2.template triangularView<Lower>(); |
| 394 | m3 = m2.template triangularView<Lower>(); |
| 395 | VERIFY_IS_APPROX(m3, refMat3); |
| 396 | |
| 397 | refMat3 = refMat2.template triangularView<Upper>(); |
| 398 | m3 = m2.template triangularView<Upper>(); |
| 399 | VERIFY_IS_APPROX(m3, refMat3); |
| 400 | |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 401 | { |
| 402 | refMat3 = refMat2.template triangularView<UnitUpper>(); |
| 403 | m3 = m2.template triangularView<UnitUpper>(); |
| 404 | VERIFY_IS_APPROX(m3, refMat3); |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 405 | |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 406 | refMat3 = refMat2.template triangularView<UnitLower>(); |
| 407 | m3 = m2.template triangularView<UnitLower>(); |
| 408 | VERIFY_IS_APPROX(m3, refMat3); |
| 409 | } |
Desire NUENTSA | 4cd8245 | 2013-06-11 14:42:29 +0200 | [diff] [blame] | 410 | |
| 411 | refMat3 = refMat2.template triangularView<StrictlyUpper>(); |
| 412 | m3 = m2.template triangularView<StrictlyUpper>(); |
| 413 | VERIFY_IS_APPROX(m3, refMat3); |
| 414 | |
| 415 | refMat3 = refMat2.template triangularView<StrictlyLower>(); |
| 416 | m3 = m2.template triangularView<StrictlyLower>(); |
| 417 | VERIFY_IS_APPROX(m3, refMat3); |
Gael Guennebaud | f6b1dee | 2015-11-04 17:02:32 +0100 | [diff] [blame] | 418 | |
| 419 | // check sparse-traingular to dense |
| 420 | refMat3 = m2.template triangularView<StrictlyUpper>(); |
| 421 | VERIFY_IS_APPROX(refMat3, DenseMatrix(refMat2.template triangularView<StrictlyUpper>())); |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 422 | } |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 423 | |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 424 | // test selfadjointView |
Gael Guennebaud | 9353bba | 2011-12-04 14:39:24 +0100 | [diff] [blame] | 425 | if(!SparseMatrixType::IsRowMajor) |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 426 | { |
| 427 | DenseMatrix refMat2(rows, rows), refMat3(rows, rows); |
| 428 | SparseMatrixType m2(rows, rows), m3(rows, rows); |
| 429 | initSparse<Scalar>(density, refMat2, m2); |
| 430 | refMat3 = refMat2.template selfadjointView<Lower>(); |
| 431 | m3 = m2.template selfadjointView<Lower>(); |
| 432 | VERIFY_IS_APPROX(m3, refMat3); |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 433 | |
Gael Guennebaud | 11b492e | 2016-12-14 17:53:47 +0100 | [diff] [blame^] | 434 | refMat3 += refMat2.template selfadjointView<Lower>(); |
| 435 | m3 += m2.template selfadjointView<Lower>(); |
| 436 | VERIFY_IS_APPROX(m3, refMat3); |
| 437 | |
| 438 | refMat3 -= refMat2.template selfadjointView<Lower>(); |
| 439 | m3 -= m2.template selfadjointView<Lower>(); |
| 440 | VERIFY_IS_APPROX(m3, refMat3); |
| 441 | |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 442 | // selfadjointView only works for square matrices: |
| 443 | SparseMatrixType m4(rows, rows+1); |
| 444 | VERIFY_RAISES_ASSERT(m4.template selfadjointView<Lower>()); |
| 445 | VERIFY_RAISES_ASSERT(m4.template selfadjointView<Upper>()); |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 446 | } |
| 447 | |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 448 | // test sparseView |
| 449 | { |
| 450 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows); |
| 451 | SparseMatrixType m2(rows, rows); |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 452 | initSparse<Scalar>(density, refMat2, m2); |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 453 | VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval()); |
Gael Guennebaud | 8456bbb | 2016-05-18 16:53:28 +0200 | [diff] [blame] | 454 | |
| 455 | // sparse view on expressions: |
| 456 | VERIFY_IS_APPROX((s1*m2).eval(), (s1*refMat2).sparseView().eval()); |
| 457 | VERIFY_IS_APPROX((m2+m2).eval(), (refMat2+refMat2).sparseView().eval()); |
| 458 | VERIFY_IS_APPROX((m2*m2).eval(), (refMat2.lazyProduct(refMat2)).sparseView().eval()); |
| 459 | VERIFY_IS_APPROX((m2*m2).eval(), (refMat2*refMat2).sparseView().eval()); |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 460 | } |
Gael Guennebaud | 82f9aa1 | 2011-12-04 21:49:21 +0100 | [diff] [blame] | 461 | |
| 462 | // test diagonal |
| 463 | { |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 464 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols); |
| 465 | SparseMatrixType m2(rows, cols); |
Gael Guennebaud | 82f9aa1 | 2011-12-04 21:49:21 +0100 | [diff] [blame] | 466 | initSparse<Scalar>(density, refMat2, m2); |
| 467 | VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval()); |
Gael Guennebaud | b26e697 | 2014-12-01 14:41:39 +0100 | [diff] [blame] | 468 | VERIFY_IS_APPROX(const_cast<const SparseMatrixType&>(m2).diagonal(), refMat2.diagonal().eval()); |
| 469 | |
| 470 | initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag); |
| 471 | m2.diagonal() += refMat2.diagonal(); |
| 472 | refMat2.diagonal() += refMat2.diagonal(); |
| 473 | VERIFY_IS_APPROX(m2, refMat2); |
Gael Guennebaud | 82f9aa1 | 2011-12-04 21:49:21 +0100 | [diff] [blame] | 474 | } |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 475 | |
Gael Guennebaud | 62f21e2 | 2015-06-24 17:55:00 +0200 | [diff] [blame] | 476 | // test diagonal to sparse |
| 477 | { |
| 478 | DenseVector d = DenseVector::Random(rows); |
| 479 | DenseMatrix refMat2 = d.asDiagonal(); |
| 480 | SparseMatrixType m2(rows, rows); |
| 481 | m2 = d.asDiagonal(); |
| 482 | VERIFY_IS_APPROX(m2, refMat2); |
Gael Guennebaud | 4c8cd13 | 2015-06-24 18:11:06 +0200 | [diff] [blame] | 483 | SparseMatrixType m3(d.asDiagonal()); |
| 484 | VERIFY_IS_APPROX(m3, refMat2); |
Gael Guennebaud | 62f21e2 | 2015-06-24 17:55:00 +0200 | [diff] [blame] | 485 | refMat2 += d.asDiagonal(); |
| 486 | m2 += d.asDiagonal(); |
| 487 | VERIFY_IS_APPROX(m2, refMat2); |
| 488 | } |
| 489 | |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 490 | // test conservative resize |
| 491 | { |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 492 | std::vector< std::pair<StorageIndex,StorageIndex> > inc; |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 493 | if(rows > 3 && cols > 2) |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 494 | inc.push_back(std::pair<StorageIndex,StorageIndex>(-3,-2)); |
| 495 | inc.push_back(std::pair<StorageIndex,StorageIndex>(0,0)); |
| 496 | inc.push_back(std::pair<StorageIndex,StorageIndex>(3,2)); |
| 497 | inc.push_back(std::pair<StorageIndex,StorageIndex>(3,0)); |
| 498 | inc.push_back(std::pair<StorageIndex,StorageIndex>(0,3)); |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 499 | |
| 500 | for(size_t i = 0; i< inc.size(); i++) { |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 501 | StorageIndex incRows = inc[i].first; |
| 502 | StorageIndex incCols = inc[i].second; |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 503 | SparseMatrixType m1(rows, cols); |
| 504 | DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols); |
| 505 | initSparse<Scalar>(density, refMat1, m1); |
| 506 | |
| 507 | m1.conservativeResize(rows+incRows, cols+incCols); |
| 508 | refMat1.conservativeResize(rows+incRows, cols+incCols); |
| 509 | if (incRows > 0) refMat1.bottomRows(incRows).setZero(); |
| 510 | if (incCols > 0) refMat1.rightCols(incCols).setZero(); |
| 511 | |
| 512 | VERIFY_IS_APPROX(m1, refMat1); |
| 513 | |
| 514 | // Insert new values |
| 515 | if (incRows > 0) |
Gael Guennebaud | 7ee378d | 2013-07-12 16:40:02 +0200 | [diff] [blame] | 516 | m1.insert(m1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1; |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 517 | if (incCols > 0) |
Gael Guennebaud | 7ee378d | 2013-07-12 16:40:02 +0200 | [diff] [blame] | 518 | m1.insert(0, m1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1; |
Benjamin Piwowarski | 6bf49ce | 2012-07-19 00:07:06 +0200 | [diff] [blame] | 519 | |
| 520 | VERIFY_IS_APPROX(m1, refMat1); |
| 521 | |
| 522 | |
| 523 | } |
| 524 | } |
Desire NUENTSA | 4cd8245 | 2013-06-11 14:42:29 +0200 | [diff] [blame] | 525 | |
| 526 | // test Identity matrix |
| 527 | { |
| 528 | DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows); |
| 529 | SparseMatrixType m1(rows, rows); |
| 530 | m1.setIdentity(); |
| 531 | VERIFY_IS_APPROX(m1, refMat1); |
Gael Guennebaud | 8a211bb | 2015-10-25 22:01:58 +0100 | [diff] [blame] | 532 | for(int k=0; k<rows*rows/4; ++k) |
| 533 | { |
| 534 | Index i = internal::random<Index>(0,rows-1); |
| 535 | Index j = internal::random<Index>(0,rows-1); |
Gael Guennebaud | 73f692d | 2015-10-27 11:01:37 +0100 | [diff] [blame] | 536 | Scalar v = internal::random<Scalar>(); |
Gael Guennebaud | 8a211bb | 2015-10-25 22:01:58 +0100 | [diff] [blame] | 537 | m1.coeffRef(i,j) = v; |
| 538 | refMat1.coeffRef(i,j) = v; |
| 539 | VERIFY_IS_APPROX(m1, refMat1); |
| 540 | if(internal::random<Index>(0,10)<2) |
| 541 | m1.makeCompressed(); |
| 542 | } |
| 543 | m1.setIdentity(); |
| 544 | refMat1.setIdentity(); |
| 545 | VERIFY_IS_APPROX(m1, refMat1); |
Desire NUENTSA | 4cd8245 | 2013-06-11 14:42:29 +0200 | [diff] [blame] | 546 | } |
Gael Guennebaud | ec46970 | 2016-02-01 15:04:33 +0100 | [diff] [blame] | 547 | |
| 548 | // test array/vector of InnerIterator |
| 549 | { |
| 550 | typedef typename SparseMatrixType::InnerIterator IteratorType; |
| 551 | |
| 552 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, cols); |
| 553 | SparseMatrixType m2(rows, cols); |
| 554 | initSparse<Scalar>(density, refMat2, m2); |
| 555 | IteratorType static_array[2]; |
| 556 | static_array[0] = IteratorType(m2,0); |
| 557 | static_array[1] = IteratorType(m2,m2.outerSize()-1); |
| 558 | VERIFY( static_array[0] || m2.innerVector(static_array[0].outer()).nonZeros() == 0 ); |
| 559 | VERIFY( static_array[1] || m2.innerVector(static_array[1].outer()).nonZeros() == 0 ); |
| 560 | if(static_array[0] && static_array[1]) |
| 561 | { |
| 562 | ++(static_array[1]); |
| 563 | static_array[1] = IteratorType(m2,0); |
| 564 | VERIFY( static_array[1] ); |
| 565 | VERIFY( static_array[1].index() == static_array[0].index() ); |
| 566 | VERIFY( static_array[1].outer() == static_array[0].outer() ); |
| 567 | VERIFY( static_array[1].value() == static_array[0].value() ); |
| 568 | } |
| 569 | |
| 570 | std::vector<IteratorType> iters(2); |
| 571 | iters[0] = IteratorType(m2,0); |
| 572 | iters[1] = IteratorType(m2,m2.outerSize()-1); |
| 573 | } |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 574 | } |
| 575 | |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 576 | |
Christoph Hertzberg | c5a3777 | 2014-10-31 17:19:05 +0100 | [diff] [blame] | 577 | template<typename SparseMatrixType> |
Christoph Hertzberg | e8cdbed | 2014-12-04 22:48:53 +0100 | [diff] [blame] | 578 | void big_sparse_triplet(Index rows, Index cols, double density) { |
| 579 | typedef typename SparseMatrixType::StorageIndex StorageIndex; |
| 580 | typedef typename SparseMatrixType::Scalar Scalar; |
| 581 | typedef Triplet<Scalar,Index> TripletType; |
| 582 | std::vector<TripletType> triplets; |
| 583 | double nelements = density * rows*cols; |
| 584 | VERIFY(nelements>=0 && nelements < NumTraits<StorageIndex>::highest()); |
| 585 | Index ntriplets = Index(nelements); |
| 586 | triplets.reserve(ntriplets); |
| 587 | Scalar sum = Scalar(0); |
| 588 | for(Index i=0;i<ntriplets;++i) |
| 589 | { |
| 590 | Index r = internal::random<Index>(0,rows-1); |
| 591 | Index c = internal::random<Index>(0,cols-1); |
| 592 | Scalar v = internal::random<Scalar>(); |
| 593 | triplets.push_back(TripletType(r,c,v)); |
| 594 | sum += v; |
| 595 | } |
| 596 | SparseMatrixType m(rows,cols); |
| 597 | m.setFromTriplets(triplets.begin(), triplets.end()); |
| 598 | VERIFY(m.nonZeros() <= ntriplets); |
| 599 | VERIFY_IS_APPROX(sum, m.sum()); |
Christoph Hertzberg | c5a3777 | 2014-10-31 17:19:05 +0100 | [diff] [blame] | 600 | } |
| 601 | |
| 602 | |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 603 | void test_sparse_basic() |
| 604 | { |
| 605 | for(int i = 0; i < g_repeat; i++) { |
Gael Guennebaud | c43154b | 2015-03-04 10:16:46 +0100 | [diff] [blame] | 606 | 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] | 607 | if(Eigen::internal::random<int>(0,4) == 0) { |
| 608 | r = c; // check square matrices in 25% of tries |
| 609 | } |
| 610 | EIGEN_UNUSED_VARIABLE(r+c); |
| 611 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(1, 1)) )); |
Gael Guennebaud | 91fe150 | 2011-06-07 11:28:16 +0200 | [diff] [blame] | 612 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) )); |
Christoph Hertzberg | 0833b82 | 2014-10-31 17:12:13 +0100 | [diff] [blame] | 613 | CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(r, c)) )); |
| 614 | CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(r, c)) )); |
| 615 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(r, c)) )); |
Gael Guennebaud | d7698c1 | 2015-03-19 15:11:05 +0100 | [diff] [blame] | 616 | CALL_SUBTEST_5(( sparse_basic(SparseMatrix<double,ColMajor,long int>(r, c)) )); |
| 617 | CALL_SUBTEST_5(( sparse_basic(SparseMatrix<double,RowMajor,long int>(r, c)) )); |
Gael Guennebaud | 47e8902 | 2013-06-10 10:34:03 +0200 | [diff] [blame] | 618 | |
Gael Guennebaud | c43154b | 2015-03-04 10:16:46 +0100 | [diff] [blame] | 619 | r = Eigen::internal::random<int>(1,100); |
| 620 | c = Eigen::internal::random<int>(1,100); |
| 621 | if(Eigen::internal::random<int>(0,4) == 0) { |
| 622 | r = c; // check square matrices in 25% of tries |
| 623 | } |
| 624 | |
Gael Guennebaud | d7698c1 | 2015-03-19 15:11:05 +0100 | [diff] [blame] | 625 | CALL_SUBTEST_6(( sparse_basic(SparseMatrix<double,ColMajor,short int>(short(r), short(c))) )); |
| 626 | 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] | 627 | } |
Christoph Hertzberg | c5a3777 | 2014-10-31 17:19:05 +0100 | [diff] [blame] | 628 | |
| 629 | // Regression test for bug 900: (manually insert higher values here, if you have enough RAM): |
| 630 | CALL_SUBTEST_3((big_sparse_triplet<SparseMatrix<float, RowMajor, int> >(10000, 10000, 0.125))); |
| 631 | 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] | 632 | |
| 633 | // Regression test for bug 1105 |
Christoph Hertzberg | 7268b10 | 2016-05-11 19:41:53 +0200 | [diff] [blame] | 634 | #ifdef EIGEN_TEST_PART_7 |
Gael Guennebaud | bfd6ee6 | 2015-11-06 15:05:37 +0100 | [diff] [blame] | 635 | { |
| 636 | int n = Eigen::internal::random<int>(200,600); |
| 637 | SparseMatrix<std::complex<double>,0, long> mat(n, n); |
| 638 | std::complex<double> val; |
| 639 | |
| 640 | for(int i=0; i<n; ++i) |
| 641 | { |
| 642 | mat.coeffRef(i, i%(n/10)) = val; |
| 643 | VERIFY(mat.data().allocatedSize()<20*n); |
| 644 | } |
| 645 | } |
| 646 | #endif |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 647 | } |