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> |
| 6 | // |
| 7 | // Eigen is free software; you can redistribute it and/or |
| 8 | // modify it under the terms of the GNU Lesser General Public |
| 9 | // License as published by the Free Software Foundation; either |
| 10 | // version 3 of the License, or (at your option) any later version. |
| 11 | // |
| 12 | // Alternatively, you can redistribute it and/or |
| 13 | // modify it under the terms of the GNU General Public License as |
| 14 | // published by the Free Software Foundation; either version 2 of |
| 15 | // the License, or (at your option) any later version. |
| 16 | // |
| 17 | // Eigen is distributed in the hope that it will be useful, but WITHOUT ANY |
| 18 | // WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS |
| 19 | // FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the |
| 20 | // GNU General Public License for more details. |
| 21 | // |
| 22 | // You should have received a copy of the GNU Lesser General Public |
| 23 | // License and a copy of the GNU General Public License along with |
| 24 | // Eigen. If not, see <http://www.gnu.org/licenses/>. |
| 25 | |
| 26 | #include "sparse.h" |
| 27 | |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 28 | template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref) |
| 29 | { |
Hauke Heibel | f1679c7 | 2010-06-20 17:37:56 +0200 | [diff] [blame] | 30 | typedef typename SparseMatrixType::Index Index; |
| 31 | |
| 32 | const Index rows = ref.rows(); |
| 33 | const Index cols = ref.cols(); |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 34 | typedef typename SparseMatrixType::Scalar Scalar; |
| 35 | enum { Flags = SparseMatrixType::Flags }; |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 36 | |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 37 | double density = std::max(8./(rows*cols), 0.01); |
| 38 | typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix; |
| 39 | typedef Matrix<Scalar,Dynamic,1> DenseVector; |
| 40 | Scalar eps = 1e-6; |
| 41 | |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 42 | SparseMatrixType m(rows, cols); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 43 | DenseMatrix refMat = DenseMatrix::Zero(rows, cols); |
| 44 | DenseVector vec1 = DenseVector::Random(rows); |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 45 | Scalar s1 = internal::random<Scalar>(); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 46 | |
| 47 | std::vector<Vector2i> zeroCoords; |
| 48 | std::vector<Vector2i> nonzeroCoords; |
| 49 | initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords); |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 50 | |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 51 | if (zeroCoords.size()==0 || nonzeroCoords.size()==0) |
| 52 | return; |
| 53 | |
| 54 | // test coeff and coeffRef |
| 55 | for (int i=0; i<(int)zeroCoords.size(); ++i) |
| 56 | { |
| 57 | VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps ); |
Hauke Heibel | 7bc8e3a | 2010-10-25 22:13:49 +0200 | [diff] [blame] | 58 | if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value) |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 59 | VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 ); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 60 | } |
| 61 | VERIFY_IS_APPROX(m, refMat); |
| 62 | |
| 63 | m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5); |
| 64 | refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5); |
| 65 | |
| 66 | VERIFY_IS_APPROX(m, refMat); |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 67 | /* |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 68 | // test InnerIterators and Block expressions |
| 69 | for (int t=0; t<10; ++t) |
| 70 | { |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 71 | int j = internal::random<int>(0,cols-1); |
| 72 | int i = internal::random<int>(0,rows-1); |
| 73 | int w = internal::random<int>(1,cols-j-1); |
| 74 | int h = internal::random<int>(1,rows-i-1); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 75 | |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 76 | // VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w)); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 77 | for(int c=0; c<w; c++) |
| 78 | { |
| 79 | VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c)); |
| 80 | for(int r=0; r<h; r++) |
| 81 | { |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 82 | // VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r)); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 83 | } |
| 84 | } |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 85 | // for(int r=0; r<h; r++) |
| 86 | // { |
| 87 | // VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r)); |
| 88 | // for(int c=0; c<w; c++) |
| 89 | // { |
| 90 | // VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c)); |
| 91 | // } |
| 92 | // } |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 93 | } |
| 94 | |
| 95 | for(int c=0; c<cols; c++) |
| 96 | { |
| 97 | VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c)); |
| 98 | VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c)); |
| 99 | } |
| 100 | |
| 101 | for(int r=0; r<rows; r++) |
| 102 | { |
| 103 | VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r)); |
| 104 | VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r)); |
| 105 | } |
| 106 | */ |
| 107 | |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 108 | // test insert (inner random) |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 109 | { |
| 110 | DenseMatrix m1(rows,cols); |
| 111 | m1.setZero(); |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 112 | SparseMatrixType m2(rows,cols); |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 113 | m2.reserve(10); |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 114 | for (int j=0; j<cols; ++j) |
| 115 | { |
| 116 | for (int k=0; k<rows/2; ++k) |
| 117 | { |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 118 | int i = internal::random<int>(0,rows-1); |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 119 | if (m1.coeff(i,j)==Scalar(0)) |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 120 | m2.insert(i,j) = m1(i,j) = internal::random<Scalar>(); |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 121 | } |
| 122 | } |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 123 | m2.finalize(); |
| 124 | VERIFY_IS_APPROX(m2,m1); |
| 125 | } |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 126 | |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 127 | // test insert (fully random) |
| 128 | { |
| 129 | DenseMatrix m1(rows,cols); |
| 130 | m1.setZero(); |
| 131 | SparseMatrixType m2(rows,cols); |
| 132 | m2.reserve(10); |
| 133 | for (int k=0; k<rows*cols; ++k) |
| 134 | { |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 135 | int i = internal::random<int>(0,rows-1); |
| 136 | int j = internal::random<int>(0,cols-1); |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 137 | if (m1.coeff(i,j)==Scalar(0)) |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 138 | m2.insert(i,j) = m1(i,j) = internal::random<Scalar>(); |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 139 | } |
| 140 | m2.finalize(); |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 141 | VERIFY_IS_APPROX(m2,m1); |
Gael Guennebaud | 5015e48 | 2008-12-11 18:26:24 +0000 | [diff] [blame] | 142 | } |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 143 | |
Gael Guennebaud | 2d53466 | 2009-01-14 21:27:54 +0000 | [diff] [blame] | 144 | // test basic computations |
| 145 | { |
| 146 | DenseMatrix refM1 = DenseMatrix::Zero(rows, rows); |
| 147 | DenseMatrix refM2 = DenseMatrix::Zero(rows, rows); |
| 148 | DenseMatrix refM3 = DenseMatrix::Zero(rows, rows); |
| 149 | DenseMatrix refM4 = DenseMatrix::Zero(rows, rows); |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 150 | SparseMatrixType m1(rows, rows); |
| 151 | SparseMatrixType m2(rows, rows); |
| 152 | SparseMatrixType m3(rows, rows); |
| 153 | SparseMatrixType m4(rows, rows); |
Gael Guennebaud | 2d53466 | 2009-01-14 21:27:54 +0000 | [diff] [blame] | 154 | initSparse<Scalar>(density, refM1, m1); |
| 155 | initSparse<Scalar>(density, refM2, m2); |
| 156 | initSparse<Scalar>(density, refM3, m3); |
| 157 | initSparse<Scalar>(density, refM4, m4); |
| 158 | |
| 159 | VERIFY_IS_APPROX(m1+m2, refM1+refM2); |
| 160 | VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3); |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 161 | VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2)); |
Gael Guennebaud | 2d53466 | 2009-01-14 21:27:54 +0000 | [diff] [blame] | 162 | VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2); |
| 163 | |
| 164 | VERIFY_IS_APPROX(m1*=s1, refM1*=s1); |
| 165 | VERIFY_IS_APPROX(m1/=s1, refM1/=s1); |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 166 | |
Gael Guennebaud | e7c48fa | 2009-01-23 13:59:32 +0000 | [diff] [blame] | 167 | VERIFY_IS_APPROX(m1+=m2, refM1+=refM2); |
| 168 | VERIFY_IS_APPROX(m1-=m2, refM1-=refM2); |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 169 | |
Gael Guennebaud | a9688f0 | 2009-02-09 09:59:30 +0000 | [diff] [blame] | 170 | VERIFY_IS_APPROX(m1.col(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0))); |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 171 | |
Gael Guennebaud | 2d53466 | 2009-01-14 21:27:54 +0000 | [diff] [blame] | 172 | refM4.setRandom(); |
| 173 | // sparse cwise* dense |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 174 | VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4)); |
Gael Guennebaud | 2d53466 | 2009-01-14 21:27:54 +0000 | [diff] [blame] | 175 | // VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4); |
| 176 | } |
| 177 | |
Gael Guennebaud | ece48a6 | 2010-06-18 11:28:30 +0200 | [diff] [blame] | 178 | // test transpose |
| 179 | { |
| 180 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows); |
| 181 | SparseMatrixType m2(rows, rows); |
| 182 | initSparse<Scalar>(density, refMat2, m2); |
| 183 | VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval()); |
| 184 | VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose()); |
| 185 | |
| 186 | VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint()); |
| 187 | } |
| 188 | |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 189 | // test innerVector() |
| 190 | { |
| 191 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows); |
Gael Guennebaud | 178858f | 2009-01-19 15:20:45 +0000 | [diff] [blame] | 192 | SparseMatrixType m2(rows, rows); |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 193 | initSparse<Scalar>(density, refMat2, m2); |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 194 | int j0 = internal::random(0,rows-1); |
| 195 | int j1 = internal::random(0,rows-1); |
Gael Guennebaud | 2d53466 | 2009-01-14 21:27:54 +0000 | [diff] [blame] | 196 | VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0)); |
| 197 | VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1)); |
Gael Guennebaud | 8ce4503 | 2009-01-27 22:48:17 +0000 | [diff] [blame] | 198 | //m2.innerVector(j0) = 2*m2.innerVector(j1); |
| 199 | //refMat2.col(j0) = 2*refMat2.col(j1); |
| 200 | //VERIFY_IS_APPROX(m2, refMat2); |
| 201 | } |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 202 | |
Gael Guennebaud | 8ce4503 | 2009-01-27 22:48:17 +0000 | [diff] [blame] | 203 | // test innerVectors() |
| 204 | { |
| 205 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows); |
| 206 | SparseMatrixType m2(rows, rows); |
| 207 | initSparse<Scalar>(density, refMat2, m2); |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 208 | int j0 = internal::random(0,rows-2); |
| 209 | int j1 = internal::random(0,rows-2); |
| 210 | int n0 = internal::random<int>(1,rows-std::max(j0,j1)); |
Gael Guennebaud | 8ce4503 | 2009-01-27 22:48:17 +0000 | [diff] [blame] | 211 | VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0)); |
| 212 | VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0), |
| 213 | refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0)); |
| 214 | //m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0); |
| 215 | //refMat2.block(0,j0,rows,n0) = refMat2.block(0,j0,rows,n0) + refMat2.block(0,j1,rows,n0); |
Gael Guennebaud | c4c7066 | 2009-01-14 14:24:10 +0000 | [diff] [blame] | 216 | } |
| 217 | |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 218 | // test prune |
| 219 | { |
| 220 | SparseMatrixType m2(rows, rows); |
| 221 | DenseMatrix refM2(rows, rows); |
| 222 | refM2.setZero(); |
| 223 | int countFalseNonZero = 0; |
| 224 | int countTrueNonZero = 0; |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 225 | for (int j=0; j<m2.outerSize(); ++j) |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 226 | { |
| 227 | m2.startVec(j); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 228 | for (int i=0; i<m2.innerSize(); ++i) |
| 229 | { |
Benoit Jacob | 4716040 | 2010-10-25 10:15:22 -0400 | [diff] [blame] | 230 | float x = internal::random<float>(0,1); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 231 | if (x<0.1) |
| 232 | { |
| 233 | // do nothing |
| 234 | } |
| 235 | else if (x<0.5) |
| 236 | { |
| 237 | countFalseNonZero++; |
Gael Guennebaud | 8710bd2 | 2010-06-02 13:32:13 +0200 | [diff] [blame] | 238 | m2.insertBackByOuterInner(j,i) = Scalar(0); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 239 | } |
| 240 | else |
| 241 | { |
| 242 | countTrueNonZero++; |
Gael Guennebaud | 8710bd2 | 2010-06-02 13:32:13 +0200 | [diff] [blame] | 243 | m2.insertBackByOuterInner(j,i) = refM2(i,j) = Scalar(1); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 244 | } |
| 245 | } |
Gael Guennebaud | 2829314 | 2009-05-04 14:25:12 +0000 | [diff] [blame] | 246 | } |
| 247 | m2.finalize(); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 248 | VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros()); |
| 249 | VERIFY_IS_APPROX(m2, refM2); |
Hauke Heibel | d204ec4 | 2010-11-02 14:33:33 +0100 | [diff] [blame] | 250 | m2.prune(Scalar(1)); |
Gael Guennebaud | 52cf07d | 2009-01-21 18:46:04 +0000 | [diff] [blame] | 251 | VERIFY(countTrueNonZero==m2.nonZeros()); |
| 252 | VERIFY_IS_APPROX(m2, refM2); |
| 253 | } |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 254 | |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 255 | // test selfadjointView |
| 256 | { |
| 257 | DenseMatrix refMat2(rows, rows), refMat3(rows, rows); |
| 258 | SparseMatrixType m2(rows, rows), m3(rows, rows); |
| 259 | initSparse<Scalar>(density, refMat2, m2); |
| 260 | refMat3 = refMat2.template selfadjointView<Lower>(); |
| 261 | m3 = m2.template selfadjointView<Lower>(); |
| 262 | VERIFY_IS_APPROX(m3, refMat3); |
| 263 | } |
| 264 | |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 265 | // test sparseView |
| 266 | { |
| 267 | DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows); |
| 268 | SparseMatrixType m2(rows, rows); |
Gael Guennebaud | 9a3ec63 | 2010-11-15 14:14:05 +0100 | [diff] [blame] | 269 | initSparse<Scalar>(density, refMat2, m2); |
Gael Guennebaud | fa6d36e | 2010-07-22 15:57:01 +0200 | [diff] [blame] | 270 | VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval()); |
| 271 | } |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 272 | } |
| 273 | |
| 274 | void test_sparse_basic() |
| 275 | { |
| 276 | for(int i = 0; i < g_repeat; i++) { |
Benoit Jacob | 2840ac7 | 2009-10-28 18:19:29 -0400 | [diff] [blame] | 277 | CALL_SUBTEST_1( sparse_basic(SparseMatrix<double>(8, 8)) ); |
| 278 | CALL_SUBTEST_2( sparse_basic(SparseMatrix<std::complex<double> >(16, 16)) ); |
| 279 | CALL_SUBTEST_1( sparse_basic(SparseMatrix<double>(33, 33)) ); |
Gael Guennebaud | 9f79558 | 2009-12-16 19:18:40 +0100 | [diff] [blame] | 280 | |
Benoit Jacob | 2840ac7 | 2009-10-28 18:19:29 -0400 | [diff] [blame] | 281 | CALL_SUBTEST_3( sparse_basic(DynamicSparseMatrix<double>(8, 8)) ); |
Gael Guennebaud | 86ccd99 | 2008-11-05 13:47:55 +0000 | [diff] [blame] | 282 | } |
| 283 | } |