| 1 | // This file is part of Eigen, a lightweight C++ template library
|
|---|
| 2 | // for linear algebra.
|
|---|
| 3 | //
|
|---|
| 4 | // Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
|
|---|
| 5 | // Copyright (C) 2008 Daniel Gomez Ferro <dgomezferro@gmail.com>
|
|---|
| 6 | // Copyright (C) 2013 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
|
|---|
| 7 | //
|
|---|
| 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/.
|
|---|
| 11 |
|
|---|
| 12 | #include "sparse.h"
|
|---|
| 13 |
|
|---|
| 14 | template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref)
|
|---|
| 15 | {
|
|---|
| 16 | typedef typename SparseMatrixType::Index Index;
|
|---|
| 17 | typedef Matrix<Index,2,1> Vector2;
|
|---|
| 18 |
|
|---|
| 19 | const Index rows = ref.rows();
|
|---|
| 20 | const Index cols = ref.cols();
|
|---|
| 21 | typedef typename SparseMatrixType::Scalar Scalar;
|
|---|
| 22 | enum { Flags = SparseMatrixType::Flags };
|
|---|
| 23 |
|
|---|
| 24 | double density = (std::max)(8./(rows*cols), 0.01);
|
|---|
| 25 | typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|---|
| 26 | typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
|---|
| 27 | typedef Matrix<Scalar,1,Dynamic> RowDenseVector;
|
|---|
| 28 | Scalar eps = 1e-6;
|
|---|
| 29 |
|
|---|
| 30 | Scalar s1 = internal::random<Scalar>();
|
|---|
| 31 | {
|
|---|
| 32 | SparseMatrixType m(rows, cols);
|
|---|
| 33 | DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
|
|---|
| 34 | DenseVector vec1 = DenseVector::Random(rows);
|
|---|
| 35 |
|
|---|
| 36 | std::vector<Vector2> zeroCoords;
|
|---|
| 37 | std::vector<Vector2> nonzeroCoords;
|
|---|
| 38 | initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
|
|---|
| 39 |
|
|---|
| 40 | if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
|
|---|
| 41 | return;
|
|---|
| 42 |
|
|---|
| 43 | // test coeff and coeffRef
|
|---|
| 44 | for (int i=0; i<(int)zeroCoords.size(); ++i)
|
|---|
| 45 | {
|
|---|
| 46 | VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
|
|---|
| 47 | if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
|
|---|
| 48 | VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
|
|---|
| 49 | }
|
|---|
| 50 | VERIFY_IS_APPROX(m, refMat);
|
|---|
| 51 |
|
|---|
| 52 | m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
|
|---|
| 53 | refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
|
|---|
| 54 |
|
|---|
| 55 | VERIFY_IS_APPROX(m, refMat);
|
|---|
| 56 |
|
|---|
| 57 | // test InnerIterators and Block expressions
|
|---|
| 58 | for (int t=0; t<10; ++t)
|
|---|
| 59 | {
|
|---|
| 60 | int j = internal::random<int>(0,cols-1);
|
|---|
| 61 | int i = internal::random<int>(0,rows-1);
|
|---|
| 62 | int w = internal::random<int>(1,cols-j-1);
|
|---|
| 63 | int h = internal::random<int>(1,rows-i-1);
|
|---|
| 64 |
|
|---|
| 65 | VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
|
|---|
| 66 | for(int c=0; c<w; c++)
|
|---|
| 67 | {
|
|---|
| 68 | VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
|
|---|
| 69 | for(int r=0; r<h; r++)
|
|---|
| 70 | {
|
|---|
| 71 | VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r));
|
|---|
| 72 | VERIFY_IS_APPROX(m.block(i,j,h,w).coeff(r,c), refMat.block(i,j,h,w).coeff(r,c));
|
|---|
| 73 | }
|
|---|
| 74 | }
|
|---|
| 75 | for(int r=0; r<h; r++)
|
|---|
| 76 | {
|
|---|
| 77 | VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
|
|---|
| 78 | for(int c=0; c<w; c++)
|
|---|
| 79 | {
|
|---|
| 80 | VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
|
|---|
| 81 | VERIFY_IS_APPROX(m.block(i,j,h,w).coeff(r,c), refMat.block(i,j,h,w).coeff(r,c));
|
|---|
| 82 | }
|
|---|
| 83 | }
|
|---|
| 84 |
|
|---|
| 85 | VERIFY_IS_APPROX(m.middleCols(j,w), refMat.middleCols(j,w));
|
|---|
| 86 | VERIFY_IS_APPROX(m.middleRows(i,h), refMat.middleRows(i,h));
|
|---|
| 87 | for(int r=0; r<h; r++)
|
|---|
| 88 | {
|
|---|
| 89 | VERIFY_IS_APPROX(m.middleCols(j,w).row(r), refMat.middleCols(j,w).row(r));
|
|---|
| 90 | VERIFY_IS_APPROX(m.middleRows(i,h).row(r), refMat.middleRows(i,h).row(r));
|
|---|
| 91 | for(int c=0; c<w; c++)
|
|---|
| 92 | {
|
|---|
| 93 | VERIFY_IS_APPROX(m.col(c).coeff(r), refMat.col(c).coeff(r));
|
|---|
| 94 | VERIFY_IS_APPROX(m.row(r).coeff(c), refMat.row(r).coeff(c));
|
|---|
| 95 |
|
|---|
| 96 | VERIFY_IS_APPROX(m.middleCols(j,w).coeff(r,c), refMat.middleCols(j,w).coeff(r,c));
|
|---|
| 97 | VERIFY_IS_APPROX(m.middleRows(i,h).coeff(r,c), refMat.middleRows(i,h).coeff(r,c));
|
|---|
| 98 | if(m.middleCols(j,w).coeff(r,c) != Scalar(0))
|
|---|
| 99 | {
|
|---|
| 100 | VERIFY_IS_APPROX(m.middleCols(j,w).coeffRef(r,c), refMat.middleCols(j,w).coeff(r,c));
|
|---|
| 101 | }
|
|---|
| 102 | if(m.middleRows(i,h).coeff(r,c) != Scalar(0))
|
|---|
| 103 | {
|
|---|
| 104 | VERIFY_IS_APPROX(m.middleRows(i,h).coeff(r,c), refMat.middleRows(i,h).coeff(r,c));
|
|---|
| 105 | }
|
|---|
| 106 | }
|
|---|
| 107 | }
|
|---|
| 108 | for(int c=0; c<w; c++)
|
|---|
| 109 | {
|
|---|
| 110 | VERIFY_IS_APPROX(m.middleCols(j,w).col(c), refMat.middleCols(j,w).col(c));
|
|---|
| 111 | VERIFY_IS_APPROX(m.middleRows(i,h).col(c), refMat.middleRows(i,h).col(c));
|
|---|
| 112 | }
|
|---|
| 113 | }
|
|---|
| 114 |
|
|---|
| 115 | for(int c=0; c<cols; c++)
|
|---|
| 116 | {
|
|---|
| 117 | VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
|
|---|
| 118 | VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
|
|---|
| 119 | }
|
|---|
| 120 |
|
|---|
| 121 | for(int r=0; r<rows; r++)
|
|---|
| 122 | {
|
|---|
| 123 | VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
|
|---|
| 124 | VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
|
|---|
| 125 | }
|
|---|
| 126 |
|
|---|
| 127 |
|
|---|
| 128 | // test assertion
|
|---|
| 129 | VERIFY_RAISES_ASSERT( m.coeffRef(-1,1) = 0 );
|
|---|
| 130 | VERIFY_RAISES_ASSERT( m.coeffRef(0,m.cols()) = 0 );
|
|---|
| 131 | }
|
|---|
| 132 |
|
|---|
| 133 | // test insert (inner random)
|
|---|
| 134 | {
|
|---|
| 135 | DenseMatrix m1(rows,cols);
|
|---|
| 136 | m1.setZero();
|
|---|
| 137 | SparseMatrixType m2(rows,cols);
|
|---|
| 138 | if(internal::random<int>()%2)
|
|---|
| 139 | m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
|
|---|
| 140 | for (Index j=0; j<cols; ++j)
|
|---|
| 141 | {
|
|---|
| 142 | for (Index k=0; k<rows/2; ++k)
|
|---|
| 143 | {
|
|---|
| 144 | Index i = internal::random<Index>(0,rows-1);
|
|---|
| 145 | if (m1.coeff(i,j)==Scalar(0))
|
|---|
| 146 | m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
|
|---|
| 147 | }
|
|---|
| 148 | }
|
|---|
| 149 | m2.finalize();
|
|---|
| 150 | VERIFY_IS_APPROX(m2,m1);
|
|---|
| 151 | }
|
|---|
| 152 |
|
|---|
| 153 | // test insert (fully random)
|
|---|
| 154 | {
|
|---|
| 155 | DenseMatrix m1(rows,cols);
|
|---|
| 156 | m1.setZero();
|
|---|
| 157 | SparseMatrixType m2(rows,cols);
|
|---|
| 158 | if(internal::random<int>()%2)
|
|---|
| 159 | m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
|
|---|
| 160 | for (int k=0; k<rows*cols; ++k)
|
|---|
| 161 | {
|
|---|
| 162 | Index i = internal::random<Index>(0,rows-1);
|
|---|
| 163 | Index j = internal::random<Index>(0,cols-1);
|
|---|
| 164 | if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
|
|---|
| 165 | m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
|
|---|
| 166 | else
|
|---|
| 167 | {
|
|---|
| 168 | Scalar v = internal::random<Scalar>();
|
|---|
| 169 | m2.coeffRef(i,j) += v;
|
|---|
| 170 | m1(i,j) += v;
|
|---|
| 171 | }
|
|---|
| 172 | }
|
|---|
| 173 | VERIFY_IS_APPROX(m2,m1);
|
|---|
| 174 | }
|
|---|
| 175 |
|
|---|
| 176 | // test insert (un-compressed)
|
|---|
| 177 | for(int mode=0;mode<4;++mode)
|
|---|
| 178 | {
|
|---|
| 179 | DenseMatrix m1(rows,cols);
|
|---|
| 180 | m1.setZero();
|
|---|
| 181 | SparseMatrixType m2(rows,cols);
|
|---|
| 182 | VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
|
|---|
| 183 | m2.reserve(r);
|
|---|
| 184 | for (int k=0; k<rows*cols; ++k)
|
|---|
| 185 | {
|
|---|
| 186 | Index i = internal::random<Index>(0,rows-1);
|
|---|
| 187 | Index j = internal::random<Index>(0,cols-1);
|
|---|
| 188 | if (m1.coeff(i,j)==Scalar(0))
|
|---|
| 189 | m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
|
|---|
| 190 | if(mode==3)
|
|---|
| 191 | m2.reserve(r);
|
|---|
| 192 | }
|
|---|
| 193 | if(internal::random<int>()%2)
|
|---|
| 194 | m2.makeCompressed();
|
|---|
| 195 | VERIFY_IS_APPROX(m2,m1);
|
|---|
| 196 | }
|
|---|
| 197 |
|
|---|
| 198 | // test innerVector()
|
|---|
| 199 | {
|
|---|
| 200 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|---|
| 201 | SparseMatrixType m2(rows, rows);
|
|---|
| 202 | initSparse<Scalar>(density, refMat2, m2);
|
|---|
| 203 | Index j0 = internal::random<Index>(0,rows-1);
|
|---|
| 204 | Index j1 = internal::random<Index>(0,rows-1);
|
|---|
| 205 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 206 | VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
|
|---|
| 207 | else
|
|---|
| 208 | VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
|
|---|
| 209 |
|
|---|
| 210 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 211 | VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
|
|---|
| 212 | else
|
|---|
| 213 | VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
|
|---|
| 214 |
|
|---|
| 215 | SparseMatrixType m3(rows,rows);
|
|---|
| 216 | m3.reserve(VectorXi::Constant(rows,rows/2));
|
|---|
| 217 | for(Index j=0; j<rows; ++j)
|
|---|
| 218 | for(Index k=0; k<j; ++k)
|
|---|
| 219 | m3.insertByOuterInner(j,k) = k+1;
|
|---|
| 220 | for(Index j=0; j<rows; ++j)
|
|---|
| 221 | {
|
|---|
| 222 | VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
|
|---|
| 223 | if(j>0)
|
|---|
| 224 | VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
|
|---|
| 225 | }
|
|---|
| 226 | m3.makeCompressed();
|
|---|
| 227 | for(Index j=0; j<rows; ++j)
|
|---|
| 228 | {
|
|---|
| 229 | VERIFY(j==numext::real(m3.innerVector(j).nonZeros()));
|
|---|
| 230 | if(j>0)
|
|---|
| 231 | VERIFY(j==numext::real(m3.innerVector(j).lastCoeff()));
|
|---|
| 232 | }
|
|---|
| 233 |
|
|---|
| 234 | //m2.innerVector(j0) = 2*m2.innerVector(j1);
|
|---|
| 235 | //refMat2.col(j0) = 2*refMat2.col(j1);
|
|---|
| 236 | //VERIFY_IS_APPROX(m2, refMat2);
|
|---|
| 237 | }
|
|---|
| 238 |
|
|---|
| 239 | // test innerVectors()
|
|---|
| 240 | {
|
|---|
| 241 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|---|
| 242 | SparseMatrixType m2(rows, rows);
|
|---|
| 243 | initSparse<Scalar>(density, refMat2, m2);
|
|---|
| 244 | if(internal::random<float>(0,1)>0.5) m2.makeCompressed();
|
|---|
| 245 |
|
|---|
| 246 | Index j0 = internal::random<Index>(0,rows-2);
|
|---|
| 247 | Index j1 = internal::random<Index>(0,rows-2);
|
|---|
| 248 | Index n0 = internal::random<Index>(1,rows-(std::max)(j0,j1));
|
|---|
| 249 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 250 | VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
|
|---|
| 251 | else
|
|---|
| 252 | VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
|
|---|
| 253 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 254 | VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
|
|---|
| 255 | refMat2.middleRows(j0,n0)+refMat2.middleRows(j1,n0));
|
|---|
| 256 | else
|
|---|
| 257 | VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
|
|---|
| 258 | refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
|
|---|
| 259 |
|
|---|
| 260 | VERIFY_IS_APPROX(m2, refMat2);
|
|---|
| 261 |
|
|---|
| 262 | m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
|
|---|
| 263 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 264 | refMat2.middleRows(j0,n0) = (refMat2.middleRows(j0,n0) + refMat2.middleRows(j1,n0)).eval();
|
|---|
| 265 | else
|
|---|
| 266 | refMat2.middleCols(j0,n0) = (refMat2.middleCols(j0,n0) + refMat2.middleCols(j1,n0)).eval();
|
|---|
| 267 |
|
|---|
| 268 | VERIFY_IS_APPROX(m2, refMat2);
|
|---|
| 269 |
|
|---|
| 270 | }
|
|---|
| 271 |
|
|---|
| 272 | // test basic computations
|
|---|
| 273 | {
|
|---|
| 274 | DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
|
|---|
| 275 | DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
|
|---|
| 276 | DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
|
|---|
| 277 | DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
|
|---|
| 278 | SparseMatrixType m1(rows, rows);
|
|---|
| 279 | SparseMatrixType m2(rows, rows);
|
|---|
| 280 | SparseMatrixType m3(rows, rows);
|
|---|
| 281 | SparseMatrixType m4(rows, rows);
|
|---|
| 282 | initSparse<Scalar>(density, refM1, m1);
|
|---|
| 283 | initSparse<Scalar>(density, refM2, m2);
|
|---|
| 284 | initSparse<Scalar>(density, refM3, m3);
|
|---|
| 285 | initSparse<Scalar>(density, refM4, m4);
|
|---|
| 286 |
|
|---|
| 287 | VERIFY_IS_APPROX(m1+m2, refM1+refM2);
|
|---|
| 288 | VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
|
|---|
| 289 | VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
|
|---|
| 290 | VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
|
|---|
| 291 |
|
|---|
| 292 | VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
|
|---|
| 293 | VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
|
|---|
| 294 |
|
|---|
| 295 | VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
|
|---|
| 296 | VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
|
|---|
| 297 |
|
|---|
| 298 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 299 | VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
|
|---|
| 300 | else
|
|---|
| 301 | VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
|
|---|
| 302 |
|
|---|
| 303 | VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
|
|---|
| 304 | VERIFY_IS_APPROX(m1.real(), refM1.real());
|
|---|
| 305 |
|
|---|
| 306 | refM4.setRandom();
|
|---|
| 307 | // sparse cwise* dense
|
|---|
| 308 | VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
|
|---|
| 309 | // dense cwise* sparse
|
|---|
| 310 | VERIFY_IS_APPROX(refM4.cwiseProduct(m3), refM4.cwiseProduct(refM3));
|
|---|
| 311 | // VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
|
|---|
| 312 |
|
|---|
| 313 | // test aliasing
|
|---|
| 314 | VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
|
|---|
| 315 | VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
|
|---|
| 316 | VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
|
|---|
| 317 | VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
|
|---|
| 318 | }
|
|---|
| 319 |
|
|---|
| 320 | // test transpose
|
|---|
| 321 | {
|
|---|
| 322 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|---|
| 323 | SparseMatrixType m2(rows, rows);
|
|---|
| 324 | initSparse<Scalar>(density, refMat2, m2);
|
|---|
| 325 | VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
|
|---|
| 326 | VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
|
|---|
| 327 |
|
|---|
| 328 | VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
|
|---|
| 329 | }
|
|---|
| 330 |
|
|---|
| 331 |
|
|---|
| 332 |
|
|---|
| 333 | // test generic blocks
|
|---|
| 334 | {
|
|---|
| 335 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|---|
| 336 | SparseMatrixType m2(rows, rows);
|
|---|
| 337 | initSparse<Scalar>(density, refMat2, m2);
|
|---|
| 338 | Index j0 = internal::random<Index>(0,rows-2);
|
|---|
| 339 | Index j1 = internal::random<Index>(0,rows-2);
|
|---|
| 340 | Index n0 = internal::random<Index>(1,rows-(std::max)(j0,j1));
|
|---|
| 341 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 342 | VERIFY_IS_APPROX(m2.block(j0,0,n0,cols), refMat2.block(j0,0,n0,cols));
|
|---|
| 343 | else
|
|---|
| 344 | VERIFY_IS_APPROX(m2.block(0,j0,rows,n0), refMat2.block(0,j0,rows,n0));
|
|---|
| 345 |
|
|---|
| 346 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 347 | VERIFY_IS_APPROX(m2.block(j0,0,n0,cols)+m2.block(j1,0,n0,cols),
|
|---|
| 348 | refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
|
|---|
| 349 | else
|
|---|
| 350 | VERIFY_IS_APPROX(m2.block(0,j0,rows,n0)+m2.block(0,j1,rows,n0),
|
|---|
| 351 | refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
|
|---|
| 352 |
|
|---|
| 353 | Index i = internal::random<Index>(0,m2.outerSize()-1);
|
|---|
| 354 | if(SparseMatrixType::IsRowMajor) {
|
|---|
| 355 | m2.innerVector(i) = m2.innerVector(i) * s1;
|
|---|
| 356 | refMat2.row(i) = refMat2.row(i) * s1;
|
|---|
| 357 | VERIFY_IS_APPROX(m2,refMat2);
|
|---|
| 358 | } else {
|
|---|
| 359 | m2.innerVector(i) = m2.innerVector(i) * s1;
|
|---|
| 360 | refMat2.col(i) = refMat2.col(i) * s1;
|
|---|
| 361 | VERIFY_IS_APPROX(m2,refMat2);
|
|---|
| 362 | }
|
|---|
| 363 |
|
|---|
| 364 | VERIFY_IS_APPROX(DenseVector(m2.col(j0)), refMat2.col(j0));
|
|---|
| 365 | VERIFY_IS_APPROX(m2.col(j0), refMat2.col(j0));
|
|---|
| 366 |
|
|---|
| 367 | VERIFY_IS_APPROX(RowDenseVector(m2.row(j0)), refMat2.row(j0));
|
|---|
| 368 | VERIFY_IS_APPROX(m2.row(j0), refMat2.row(j0));
|
|---|
| 369 |
|
|---|
| 370 | VERIFY_IS_APPROX(m2.block(j0,j1,n0,n0), refMat2.block(j0,j1,n0,n0));
|
|---|
| 371 | VERIFY_IS_APPROX((2*m2).block(j0,j1,n0,n0), (2*refMat2).block(j0,j1,n0,n0));
|
|---|
| 372 | }
|
|---|
| 373 |
|
|---|
| 374 | // test prune
|
|---|
| 375 | {
|
|---|
| 376 | SparseMatrixType m2(rows, rows);
|
|---|
| 377 | DenseMatrix refM2(rows, rows);
|
|---|
| 378 | refM2.setZero();
|
|---|
| 379 | int countFalseNonZero = 0;
|
|---|
| 380 | int countTrueNonZero = 0;
|
|---|
| 381 | for (Index j=0; j<m2.outerSize(); ++j)
|
|---|
| 382 | {
|
|---|
| 383 | m2.startVec(j);
|
|---|
| 384 | for (Index i=0; i<m2.innerSize(); ++i)
|
|---|
| 385 | {
|
|---|
| 386 | float x = internal::random<float>(0,1);
|
|---|
| 387 | if (x<0.1)
|
|---|
| 388 | {
|
|---|
| 389 | // do nothing
|
|---|
| 390 | }
|
|---|
| 391 | else if (x<0.5)
|
|---|
| 392 | {
|
|---|
| 393 | countFalseNonZero++;
|
|---|
| 394 | m2.insertBackByOuterInner(j,i) = Scalar(0);
|
|---|
| 395 | }
|
|---|
| 396 | else
|
|---|
| 397 | {
|
|---|
| 398 | countTrueNonZero++;
|
|---|
| 399 | m2.insertBackByOuterInner(j,i) = Scalar(1);
|
|---|
| 400 | if(SparseMatrixType::IsRowMajor)
|
|---|
| 401 | refM2(j,i) = Scalar(1);
|
|---|
| 402 | else
|
|---|
| 403 | refM2(i,j) = Scalar(1);
|
|---|
| 404 | }
|
|---|
| 405 | }
|
|---|
| 406 | }
|
|---|
| 407 | m2.finalize();
|
|---|
| 408 | VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
|
|---|
| 409 | VERIFY_IS_APPROX(m2, refM2);
|
|---|
| 410 | m2.prune(Scalar(1));
|
|---|
| 411 | VERIFY(countTrueNonZero==m2.nonZeros());
|
|---|
| 412 | VERIFY_IS_APPROX(m2, refM2);
|
|---|
| 413 | }
|
|---|
| 414 |
|
|---|
| 415 | // test setFromTriplets
|
|---|
| 416 | {
|
|---|
| 417 | typedef Triplet<Scalar,Index> TripletType;
|
|---|
| 418 | std::vector<TripletType> triplets;
|
|---|
| 419 | int ntriplets = rows*cols;
|
|---|
| 420 | triplets.reserve(ntriplets);
|
|---|
| 421 | DenseMatrix refMat(rows,cols);
|
|---|
| 422 | refMat.setZero();
|
|---|
| 423 | for(int i=0;i<ntriplets;++i)
|
|---|
| 424 | {
|
|---|
| 425 | Index r = internal::random<Index>(0,rows-1);
|
|---|
| 426 | Index c = internal::random<Index>(0,cols-1);
|
|---|
| 427 | Scalar v = internal::random<Scalar>();
|
|---|
| 428 | triplets.push_back(TripletType(r,c,v));
|
|---|
| 429 | refMat(r,c) += v;
|
|---|
| 430 | }
|
|---|
| 431 | SparseMatrixType m(rows,cols);
|
|---|
| 432 | m.setFromTriplets(triplets.begin(), triplets.end());
|
|---|
| 433 | VERIFY_IS_APPROX(m, refMat);
|
|---|
| 434 | }
|
|---|
| 435 |
|
|---|
| 436 | // test triangularView
|
|---|
| 437 | {
|
|---|
| 438 | DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
|
|---|
| 439 | SparseMatrixType m2(rows, rows), m3(rows, rows);
|
|---|
| 440 | initSparse<Scalar>(density, refMat2, m2);
|
|---|
| 441 | refMat3 = refMat2.template triangularView<Lower>();
|
|---|
| 442 | m3 = m2.template triangularView<Lower>();
|
|---|
| 443 | VERIFY_IS_APPROX(m3, refMat3);
|
|---|
| 444 |
|
|---|
| 445 | refMat3 = refMat2.template triangularView<Upper>();
|
|---|
| 446 | m3 = m2.template triangularView<Upper>();
|
|---|
| 447 | VERIFY_IS_APPROX(m3, refMat3);
|
|---|
| 448 |
|
|---|
| 449 | refMat3 = refMat2.template triangularView<UnitUpper>();
|
|---|
| 450 | m3 = m2.template triangularView<UnitUpper>();
|
|---|
| 451 | VERIFY_IS_APPROX(m3, refMat3);
|
|---|
| 452 |
|
|---|
| 453 | refMat3 = refMat2.template triangularView<UnitLower>();
|
|---|
| 454 | m3 = m2.template triangularView<UnitLower>();
|
|---|
| 455 | VERIFY_IS_APPROX(m3, refMat3);
|
|---|
| 456 |
|
|---|
| 457 | refMat3 = refMat2.template triangularView<StrictlyUpper>();
|
|---|
| 458 | m3 = m2.template triangularView<StrictlyUpper>();
|
|---|
| 459 | VERIFY_IS_APPROX(m3, refMat3);
|
|---|
| 460 |
|
|---|
| 461 | refMat3 = refMat2.template triangularView<StrictlyLower>();
|
|---|
| 462 | m3 = m2.template triangularView<StrictlyLower>();
|
|---|
| 463 | VERIFY_IS_APPROX(m3, refMat3);
|
|---|
| 464 | }
|
|---|
| 465 |
|
|---|
| 466 | // test selfadjointView
|
|---|
| 467 | if(!SparseMatrixType::IsRowMajor)
|
|---|
| 468 | {
|
|---|
| 469 | DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
|
|---|
| 470 | SparseMatrixType m2(rows, rows), m3(rows, rows);
|
|---|
| 471 | initSparse<Scalar>(density, refMat2, m2);
|
|---|
| 472 | refMat3 = refMat2.template selfadjointView<Lower>();
|
|---|
| 473 | m3 = m2.template selfadjointView<Lower>();
|
|---|
| 474 | VERIFY_IS_APPROX(m3, refMat3);
|
|---|
| 475 | }
|
|---|
| 476 |
|
|---|
| 477 | // test sparseView
|
|---|
| 478 | {
|
|---|
| 479 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|---|
| 480 | SparseMatrixType m2(rows, rows);
|
|---|
| 481 | initSparse<Scalar>(density, refMat2, m2);
|
|---|
| 482 | VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
|
|---|
| 483 | }
|
|---|
| 484 |
|
|---|
| 485 | // test diagonal
|
|---|
| 486 | {
|
|---|
| 487 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|---|
| 488 | SparseMatrixType m2(rows, rows);
|
|---|
| 489 | initSparse<Scalar>(density, refMat2, m2);
|
|---|
| 490 | VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
|
|---|
| 491 | }
|
|---|
| 492 |
|
|---|
| 493 | // test conservative resize
|
|---|
| 494 | {
|
|---|
| 495 | std::vector< std::pair<Index,Index> > inc;
|
|---|
| 496 | inc.push_back(std::pair<Index,Index>(-3,-2));
|
|---|
| 497 | inc.push_back(std::pair<Index,Index>(0,0));
|
|---|
| 498 | inc.push_back(std::pair<Index,Index>(3,2));
|
|---|
| 499 | inc.push_back(std::pair<Index,Index>(3,0));
|
|---|
| 500 | inc.push_back(std::pair<Index,Index>(0,3));
|
|---|
| 501 |
|
|---|
| 502 | for(size_t i = 0; i< inc.size(); i++) {
|
|---|
| 503 | Index incRows = inc[i].first;
|
|---|
| 504 | Index incCols = inc[i].second;
|
|---|
| 505 | SparseMatrixType m1(rows, cols);
|
|---|
| 506 | DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
|
|---|
| 507 | initSparse<Scalar>(density, refMat1, m1);
|
|---|
| 508 |
|
|---|
| 509 | m1.conservativeResize(rows+incRows, cols+incCols);
|
|---|
| 510 | refMat1.conservativeResize(rows+incRows, cols+incCols);
|
|---|
| 511 | if (incRows > 0) refMat1.bottomRows(incRows).setZero();
|
|---|
| 512 | if (incCols > 0) refMat1.rightCols(incCols).setZero();
|
|---|
| 513 |
|
|---|
| 514 | VERIFY_IS_APPROX(m1, refMat1);
|
|---|
| 515 |
|
|---|
| 516 | // Insert new values
|
|---|
| 517 | if (incRows > 0)
|
|---|
| 518 | m1.insert(m1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
|
|---|
| 519 | if (incCols > 0)
|
|---|
| 520 | m1.insert(0, m1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
|
|---|
| 521 |
|
|---|
| 522 | VERIFY_IS_APPROX(m1, refMat1);
|
|---|
| 523 |
|
|---|
| 524 |
|
|---|
| 525 | }
|
|---|
| 526 | }
|
|---|
| 527 |
|
|---|
| 528 | // test Identity matrix
|
|---|
| 529 | {
|
|---|
| 530 | DenseMatrix refMat1 = DenseMatrix::Identity(rows, rows);
|
|---|
| 531 | SparseMatrixType m1(rows, rows);
|
|---|
| 532 | m1.setIdentity();
|
|---|
| 533 | VERIFY_IS_APPROX(m1, refMat1);
|
|---|
| 534 | for(int k=0; k<rows*rows/4; ++k)
|
|---|
| 535 | {
|
|---|
| 536 | Index i = internal::random<Index>(0,rows-1);
|
|---|
| 537 | Index j = internal::random<Index>(0,rows-1);
|
|---|
| 538 | Scalar v = internal::random<Scalar>();
|
|---|
| 539 | m1.coeffRef(i,j) = v;
|
|---|
| 540 | refMat1.coeffRef(i,j) = v;
|
|---|
| 541 | VERIFY_IS_APPROX(m1, refMat1);
|
|---|
| 542 | if(internal::random<Index>(0,10)<2)
|
|---|
| 543 | m1.makeCompressed();
|
|---|
| 544 | }
|
|---|
| 545 | m1.setIdentity();
|
|---|
| 546 | refMat1.setIdentity();
|
|---|
| 547 | VERIFY_IS_APPROX(m1, refMat1);
|
|---|
| 548 | }
|
|---|
| 549 | }
|
|---|
| 550 |
|
|---|
| 551 | void test_sparse_basic()
|
|---|
| 552 | {
|
|---|
| 553 | for(int i = 0; i < g_repeat; i++) {
|
|---|
| 554 | int s = Eigen::internal::random<int>(1,50);
|
|---|
| 555 | EIGEN_UNUSED_VARIABLE(s);
|
|---|
| 556 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
|
|---|
| 557 | CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(s, s)) ));
|
|---|
| 558 | CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(s, s)) ));
|
|---|
| 559 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
|
|---|
| 560 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
|
|---|
| 561 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,long int>(s, s)) ));
|
|---|
| 562 |
|
|---|
| 563 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,short int>(short(s), short(s))) ));
|
|---|
| 564 | CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,short int>(short(s), short(s))) ));
|
|---|
| 565 | }
|
|---|
| 566 | }
|
|---|