2008-08-22 01:02:47 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
2009-05-23 02:25:33 +08:00
|
|
|
// for linear algebra.
|
2008-08-22 01:02:47 +08:00
|
|
|
//
|
2011-03-22 18:58:22 +08:00
|
|
|
// Copyright (C) 2008-2011 Gael Guennebaud <gael.guennebaud@inria.fr>
|
2008-08-22 01:02:47 +08:00
|
|
|
// Copyright (C) 2008 Daniel Gomez Ferro <dgomezferro@gmail.com>
|
|
|
|
//
|
2012-07-14 02:42:47 +08:00
|
|
|
// This Source Code Form is subject to the terms of the Mozilla
|
|
|
|
// Public License v. 2.0. If a copy of the MPL was not distributed
|
|
|
|
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
|
2008-08-22 01:02:47 +08:00
|
|
|
|
2008-11-05 21:47:55 +08:00
|
|
|
#include "sparse.h"
|
2008-09-03 03:55:26 +08:00
|
|
|
|
2009-01-19 23:20:45 +08:00
|
|
|
template<typename SparseMatrixType> void sparse_basic(const SparseMatrixType& ref)
|
|
|
|
{
|
2010-06-20 23:37:56 +08:00
|
|
|
typedef typename SparseMatrixType::Index Index;
|
|
|
|
|
|
|
|
const Index rows = ref.rows();
|
|
|
|
const Index cols = ref.cols();
|
2009-01-19 23:20:45 +08:00
|
|
|
typedef typename SparseMatrixType::Scalar Scalar;
|
|
|
|
enum { Flags = SparseMatrixType::Flags };
|
2009-12-17 02:18:40 +08:00
|
|
|
|
2011-08-19 20:18:05 +08:00
|
|
|
double density = (std::max)(8./(rows*cols), 0.01);
|
2008-09-03 03:55:26 +08:00
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
|
|
|
Scalar eps = 1e-6;
|
|
|
|
|
2009-01-19 23:20:45 +08:00
|
|
|
SparseMatrixType m(rows, cols);
|
2008-09-03 03:55:26 +08:00
|
|
|
DenseMatrix refMat = DenseMatrix::Zero(rows, cols);
|
|
|
|
DenseVector vec1 = DenseVector::Random(rows);
|
2010-10-25 22:15:22 +08:00
|
|
|
Scalar s1 = internal::random<Scalar>();
|
2008-09-03 03:55:26 +08:00
|
|
|
|
|
|
|
std::vector<Vector2i> zeroCoords;
|
|
|
|
std::vector<Vector2i> nonzeroCoords;
|
|
|
|
initSparse<Scalar>(density, refMat, m, 0, &zeroCoords, &nonzeroCoords);
|
2009-12-17 02:18:40 +08:00
|
|
|
|
2008-11-05 21:47:55 +08:00
|
|
|
if (zeroCoords.size()==0 || nonzeroCoords.size()==0)
|
|
|
|
return;
|
2008-08-22 02:40:56 +08:00
|
|
|
|
|
|
|
// test coeff and coeffRef
|
2008-08-22 09:19:53 +08:00
|
|
|
for (int i=0; i<(int)zeroCoords.size(); ++i)
|
2008-08-22 02:40:56 +08:00
|
|
|
{
|
|
|
|
VERIFY_IS_MUCH_SMALLER_THAN( m.coeff(zeroCoords[i].x(),zeroCoords[i].y()), eps );
|
2010-10-26 04:13:49 +08:00
|
|
|
if(internal::is_same<SparseMatrixType,SparseMatrix<Scalar,Flags> >::value)
|
2009-01-19 23:20:45 +08:00
|
|
|
VERIFY_RAISES_ASSERT( m.coeffRef(zeroCoords[0].x(),zeroCoords[0].y()) = 5 );
|
2008-08-22 02:40:56 +08:00
|
|
|
}
|
|
|
|
VERIFY_IS_APPROX(m, refMat);
|
|
|
|
|
|
|
|
m.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
|
|
|
|
refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(m, refMat);
|
2009-01-14 22:24:10 +08:00
|
|
|
/*
|
2008-09-02 23:28:49 +08:00
|
|
|
// test InnerIterators and Block expressions
|
2008-10-21 01:03:09 +08:00
|
|
|
for (int t=0; t<10; ++t)
|
2008-09-02 23:28:49 +08:00
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
int j = internal::random<int>(0,cols-1);
|
|
|
|
int i = internal::random<int>(0,rows-1);
|
|
|
|
int w = internal::random<int>(1,cols-j-1);
|
|
|
|
int h = internal::random<int>(1,rows-i-1);
|
2008-10-21 01:03:09 +08:00
|
|
|
|
2009-01-14 22:24:10 +08:00
|
|
|
// VERIFY_IS_APPROX(m.block(i,j,h,w), refMat.block(i,j,h,w));
|
2008-10-21 01:03:09 +08:00
|
|
|
for(int c=0; c<w; c++)
|
2008-09-02 23:28:49 +08:00
|
|
|
{
|
2008-10-21 01:03:09 +08:00
|
|
|
VERIFY_IS_APPROX(m.block(i,j,h,w).col(c), refMat.block(i,j,h,w).col(c));
|
|
|
|
for(int r=0; r<h; r++)
|
2008-09-02 23:28:49 +08:00
|
|
|
{
|
2009-01-14 22:24:10 +08:00
|
|
|
// VERIFY_IS_APPROX(m.block(i,j,h,w).col(c).coeff(r), refMat.block(i,j,h,w).col(c).coeff(r));
|
2008-09-02 23:28:49 +08:00
|
|
|
}
|
|
|
|
}
|
2009-01-14 22:24:10 +08:00
|
|
|
// for(int r=0; r<h; r++)
|
|
|
|
// {
|
|
|
|
// VERIFY_IS_APPROX(m.block(i,j,h,w).row(r), refMat.block(i,j,h,w).row(r));
|
|
|
|
// for(int c=0; c<w; c++)
|
|
|
|
// {
|
|
|
|
// VERIFY_IS_APPROX(m.block(i,j,h,w).row(r).coeff(c), refMat.block(i,j,h,w).row(r).coeff(c));
|
|
|
|
// }
|
|
|
|
// }
|
2008-09-02 23:28:49 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
for(int c=0; c<cols; c++)
|
|
|
|
{
|
|
|
|
VERIFY_IS_APPROX(m.col(c) + m.col(c), (m + m).col(c));
|
|
|
|
VERIFY_IS_APPROX(m.col(c) + m.col(c), refMat.col(c) + refMat.col(c));
|
|
|
|
}
|
|
|
|
|
|
|
|
for(int r=0; r<rows; r++)
|
|
|
|
{
|
|
|
|
VERIFY_IS_APPROX(m.row(r) + m.row(r), (m + m).row(r));
|
|
|
|
VERIFY_IS_APPROX(m.row(r) + m.row(r), refMat.row(r) + refMat.row(r));
|
|
|
|
}
|
2008-11-05 21:47:55 +08:00
|
|
|
*/
|
2008-09-02 23:28:49 +08:00
|
|
|
|
2009-05-04 22:25:12 +08:00
|
|
|
// test insert (inner random)
|
2008-12-12 02:26:24 +08:00
|
|
|
{
|
|
|
|
DenseMatrix m1(rows,cols);
|
|
|
|
m1.setZero();
|
2009-01-19 23:20:45 +08:00
|
|
|
SparseMatrixType m2(rows,cols);
|
2011-12-03 02:00:16 +08:00
|
|
|
if(internal::random<int>()%2)
|
|
|
|
m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
|
2008-12-12 02:26:24 +08:00
|
|
|
for (int j=0; j<cols; ++j)
|
|
|
|
{
|
|
|
|
for (int k=0; k<rows/2; ++k)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
int i = internal::random<int>(0,rows-1);
|
2008-12-12 02:26:24 +08:00
|
|
|
if (m1.coeff(i,j)==Scalar(0))
|
2010-10-25 22:15:22 +08:00
|
|
|
m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
|
2008-12-12 02:26:24 +08:00
|
|
|
}
|
|
|
|
}
|
2009-05-04 22:25:12 +08:00
|
|
|
m2.finalize();
|
|
|
|
VERIFY_IS_APPROX(m2,m1);
|
|
|
|
}
|
2009-12-17 02:18:40 +08:00
|
|
|
|
2009-05-04 22:25:12 +08:00
|
|
|
// test insert (fully random)
|
|
|
|
{
|
|
|
|
DenseMatrix m1(rows,cols);
|
|
|
|
m1.setZero();
|
|
|
|
SparseMatrixType m2(rows,cols);
|
2011-12-03 02:00:16 +08:00
|
|
|
if(internal::random<int>()%2)
|
|
|
|
m2.reserve(VectorXi::Constant(m2.outerSize(), 2));
|
2009-05-04 22:25:12 +08:00
|
|
|
for (int k=0; k<rows*cols; ++k)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
int i = internal::random<int>(0,rows-1);
|
|
|
|
int j = internal::random<int>(0,cols-1);
|
2011-12-03 02:00:16 +08:00
|
|
|
if ((m1.coeff(i,j)==Scalar(0)) && (internal::random<int>()%2))
|
2010-10-25 22:15:22 +08:00
|
|
|
m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
|
2011-12-03 02:00:16 +08:00
|
|
|
else
|
|
|
|
{
|
|
|
|
Scalar v = internal::random<Scalar>();
|
|
|
|
m2.coeffRef(i,j) += v;
|
|
|
|
m1(i,j) += v;
|
|
|
|
}
|
2009-05-04 22:25:12 +08:00
|
|
|
}
|
2009-01-14 22:24:10 +08:00
|
|
|
VERIFY_IS_APPROX(m2,m1);
|
2008-12-12 02:26:24 +08:00
|
|
|
}
|
2011-09-08 19:42:54 +08:00
|
|
|
|
|
|
|
// test insert (un-compressed)
|
|
|
|
for(int mode=0;mode<4;++mode)
|
|
|
|
{
|
|
|
|
DenseMatrix m1(rows,cols);
|
|
|
|
m1.setZero();
|
|
|
|
SparseMatrixType m2(rows,cols);
|
|
|
|
VectorXi r(VectorXi::Constant(m2.outerSize(), ((mode%2)==0) ? m2.innerSize() : std::max<int>(1,m2.innerSize()/8)));
|
|
|
|
m2.reserve(r);
|
|
|
|
for (int k=0; k<rows*cols; ++k)
|
|
|
|
{
|
|
|
|
int i = internal::random<int>(0,rows-1);
|
|
|
|
int j = internal::random<int>(0,cols-1);
|
|
|
|
if (m1.coeff(i,j)==Scalar(0))
|
|
|
|
m2.insert(i,j) = m1(i,j) = internal::random<Scalar>();
|
|
|
|
if(mode==3)
|
|
|
|
m2.reserve(r);
|
|
|
|
}
|
2011-12-03 02:00:16 +08:00
|
|
|
if(internal::random<int>()%2)
|
|
|
|
m2.makeCompressed();
|
2011-09-08 19:42:54 +08:00
|
|
|
VERIFY_IS_APPROX(m2,m1);
|
|
|
|
}
|
2009-12-17 02:18:40 +08:00
|
|
|
|
2009-01-15 05:27:54 +08:00
|
|
|
// test basic computations
|
|
|
|
{
|
|
|
|
DenseMatrix refM1 = DenseMatrix::Zero(rows, rows);
|
|
|
|
DenseMatrix refM2 = DenseMatrix::Zero(rows, rows);
|
|
|
|
DenseMatrix refM3 = DenseMatrix::Zero(rows, rows);
|
|
|
|
DenseMatrix refM4 = DenseMatrix::Zero(rows, rows);
|
2009-01-19 23:20:45 +08:00
|
|
|
SparseMatrixType m1(rows, rows);
|
|
|
|
SparseMatrixType m2(rows, rows);
|
|
|
|
SparseMatrixType m3(rows, rows);
|
|
|
|
SparseMatrixType m4(rows, rows);
|
2009-01-15 05:27:54 +08:00
|
|
|
initSparse<Scalar>(density, refM1, m1);
|
|
|
|
initSparse<Scalar>(density, refM2, m2);
|
|
|
|
initSparse<Scalar>(density, refM3, m3);
|
|
|
|
initSparse<Scalar>(density, refM4, m4);
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(m1+m2, refM1+refM2);
|
|
|
|
VERIFY_IS_APPROX(m1+m2+m3, refM1+refM2+refM3);
|
2009-12-17 02:18:40 +08:00
|
|
|
VERIFY_IS_APPROX(m3.cwiseProduct(m1+m2), refM3.cwiseProduct(refM1+refM2));
|
2009-01-15 05:27:54 +08:00
|
|
|
VERIFY_IS_APPROX(m1*s1-m2, refM1*s1-refM2);
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(m1*=s1, refM1*=s1);
|
|
|
|
VERIFY_IS_APPROX(m1/=s1, refM1/=s1);
|
2009-12-17 02:18:40 +08:00
|
|
|
|
2009-01-23 21:59:32 +08:00
|
|
|
VERIFY_IS_APPROX(m1+=m2, refM1+=refM2);
|
|
|
|
VERIFY_IS_APPROX(m1-=m2, refM1-=refM2);
|
2009-12-17 02:18:40 +08:00
|
|
|
|
2011-12-04 21:39:24 +08:00
|
|
|
if(SparseMatrixType::IsRowMajor)
|
|
|
|
VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.row(0).dot(refM2.row(0)));
|
|
|
|
else
|
|
|
|
VERIFY_IS_APPROX(m1.innerVector(0).dot(refM2.row(0)), refM1.col(0).dot(refM2.row(0)));
|
2009-12-17 02:18:40 +08:00
|
|
|
|
2011-12-04 06:49:37 +08:00
|
|
|
VERIFY_IS_APPROX(m1.conjugate(), refM1.conjugate());
|
|
|
|
VERIFY_IS_APPROX(m1.real(), refM1.real());
|
|
|
|
|
2009-01-15 05:27:54 +08:00
|
|
|
refM4.setRandom();
|
|
|
|
// sparse cwise* dense
|
2009-12-17 02:18:40 +08:00
|
|
|
VERIFY_IS_APPROX(m3.cwiseProduct(refM4), refM3.cwiseProduct(refM4));
|
2009-01-15 05:27:54 +08:00
|
|
|
// VERIFY_IS_APPROX(m3.cwise()/refM4, refM3.cwise()/refM4);
|
2012-07-25 15:33:50 +08:00
|
|
|
|
|
|
|
// test aliasing
|
|
|
|
VERIFY_IS_APPROX((m1 = -m1), (refM1 = -refM1));
|
|
|
|
VERIFY_IS_APPROX((m1 = m1.transpose()), (refM1 = refM1.transpose().eval()));
|
|
|
|
VERIFY_IS_APPROX((m1 = -m1.transpose()), (refM1 = -refM1.transpose().eval()));
|
|
|
|
VERIFY_IS_APPROX((m1 += -m1), (refM1 += -refM1));
|
2009-01-15 05:27:54 +08:00
|
|
|
}
|
|
|
|
|
2010-06-18 17:28:30 +08:00
|
|
|
// test transpose
|
|
|
|
{
|
|
|
|
DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|
|
|
SparseMatrixType m2(rows, rows);
|
|
|
|
initSparse<Scalar>(density, refMat2, m2);
|
|
|
|
VERIFY_IS_APPROX(m2.transpose().eval(), refMat2.transpose().eval());
|
|
|
|
VERIFY_IS_APPROX(m2.transpose(), refMat2.transpose());
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(SparseMatrixType(m2.adjoint()), refMat2.adjoint());
|
|
|
|
}
|
|
|
|
|
2009-01-14 22:24:10 +08:00
|
|
|
// test innerVector()
|
|
|
|
{
|
|
|
|
DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
2009-01-19 23:20:45 +08:00
|
|
|
SparseMatrixType m2(rows, rows);
|
2009-01-14 22:24:10 +08:00
|
|
|
initSparse<Scalar>(density, refMat2, m2);
|
2011-06-07 17:28:16 +08:00
|
|
|
int j0 = internal::random<int>(0,rows-1);
|
|
|
|
int j1 = internal::random<int>(0,rows-1);
|
2011-12-04 21:39:24 +08:00
|
|
|
if(SparseMatrixType::IsRowMajor)
|
|
|
|
VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.row(j0));
|
|
|
|
else
|
|
|
|
VERIFY_IS_APPROX(m2.innerVector(j0), refMat2.col(j0));
|
|
|
|
|
|
|
|
if(SparseMatrixType::IsRowMajor)
|
|
|
|
VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.row(j0)+refMat2.row(j1));
|
|
|
|
else
|
|
|
|
VERIFY_IS_APPROX(m2.innerVector(j0)+m2.innerVector(j1), refMat2.col(j0)+refMat2.col(j1));
|
2011-12-21 01:10:02 +08:00
|
|
|
|
|
|
|
SparseMatrixType m3(rows,rows);
|
|
|
|
m3.reserve(VectorXi::Constant(rows,rows/2));
|
|
|
|
for(int j=0; j<rows; ++j)
|
|
|
|
for(int k=0; k<j; ++k)
|
|
|
|
m3.insertByOuterInner(j,k) = k+1;
|
|
|
|
for(int j=0; j<rows; ++j)
|
|
|
|
{
|
2011-12-23 16:41:14 +08:00
|
|
|
VERIFY(j==internal::real(m3.innerVector(j).nonZeros()));
|
2011-12-21 01:10:02 +08:00
|
|
|
if(j>0)
|
2011-12-23 16:41:14 +08:00
|
|
|
VERIFY(j==internal::real(m3.innerVector(j).lastCoeff()));
|
2011-12-21 01:10:02 +08:00
|
|
|
}
|
|
|
|
m3.makeCompressed();
|
|
|
|
for(int j=0; j<rows; ++j)
|
|
|
|
{
|
2011-12-23 16:41:14 +08:00
|
|
|
VERIFY(j==internal::real(m3.innerVector(j).nonZeros()));
|
2011-12-21 01:10:02 +08:00
|
|
|
if(j>0)
|
2011-12-23 16:41:14 +08:00
|
|
|
VERIFY(j==internal::real(m3.innerVector(j).lastCoeff()));
|
2011-12-21 01:10:02 +08:00
|
|
|
}
|
|
|
|
|
2009-01-28 06:48:17 +08:00
|
|
|
//m2.innerVector(j0) = 2*m2.innerVector(j1);
|
|
|
|
//refMat2.col(j0) = 2*refMat2.col(j1);
|
|
|
|
//VERIFY_IS_APPROX(m2, refMat2);
|
|
|
|
}
|
2009-12-17 02:18:40 +08:00
|
|
|
|
2009-01-28 06:48:17 +08:00
|
|
|
// test innerVectors()
|
|
|
|
{
|
|
|
|
DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|
|
|
SparseMatrixType m2(rows, rows);
|
|
|
|
initSparse<Scalar>(density, refMat2, m2);
|
2011-06-07 17:28:16 +08:00
|
|
|
int j0 = internal::random<int>(0,rows-2);
|
|
|
|
int j1 = internal::random<int>(0,rows-2);
|
2011-08-19 20:18:05 +08:00
|
|
|
int n0 = internal::random<int>(1,rows-(std::max)(j0,j1));
|
2011-12-04 21:39:24 +08:00
|
|
|
if(SparseMatrixType::IsRowMajor)
|
|
|
|
VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(j0,0,n0,cols));
|
|
|
|
else
|
|
|
|
VERIFY_IS_APPROX(m2.innerVectors(j0,n0), refMat2.block(0,j0,rows,n0));
|
|
|
|
if(SparseMatrixType::IsRowMajor)
|
|
|
|
VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
|
|
|
|
refMat2.block(j0,0,n0,cols)+refMat2.block(j1,0,n0,cols));
|
|
|
|
else
|
|
|
|
VERIFY_IS_APPROX(m2.innerVectors(j0,n0)+m2.innerVectors(j1,n0),
|
|
|
|
refMat2.block(0,j0,rows,n0)+refMat2.block(0,j1,rows,n0));
|
2009-01-28 06:48:17 +08:00
|
|
|
//m2.innerVectors(j0,n0) = m2.innerVectors(j0,n0) + m2.innerVectors(j1,n0);
|
|
|
|
//refMat2.block(0,j0,rows,n0) = refMat2.block(0,j0,rows,n0) + refMat2.block(0,j1,rows,n0);
|
2009-01-14 22:24:10 +08:00
|
|
|
}
|
|
|
|
|
2009-01-22 02:46:04 +08:00
|
|
|
// test prune
|
|
|
|
{
|
|
|
|
SparseMatrixType m2(rows, rows);
|
|
|
|
DenseMatrix refM2(rows, rows);
|
|
|
|
refM2.setZero();
|
|
|
|
int countFalseNonZero = 0;
|
|
|
|
int countTrueNonZero = 0;
|
|
|
|
for (int j=0; j<m2.outerSize(); ++j)
|
2009-05-04 22:25:12 +08:00
|
|
|
{
|
|
|
|
m2.startVec(j);
|
2009-01-22 02:46:04 +08:00
|
|
|
for (int i=0; i<m2.innerSize(); ++i)
|
|
|
|
{
|
2010-10-25 22:15:22 +08:00
|
|
|
float x = internal::random<float>(0,1);
|
2009-01-22 02:46:04 +08:00
|
|
|
if (x<0.1)
|
|
|
|
{
|
|
|
|
// do nothing
|
|
|
|
}
|
|
|
|
else if (x<0.5)
|
|
|
|
{
|
|
|
|
countFalseNonZero++;
|
2010-06-02 19:32:13 +08:00
|
|
|
m2.insertBackByOuterInner(j,i) = Scalar(0);
|
2009-01-22 02:46:04 +08:00
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
countTrueNonZero++;
|
2011-12-04 21:39:24 +08:00
|
|
|
m2.insertBackByOuterInner(j,i) = Scalar(1);
|
|
|
|
if(SparseMatrixType::IsRowMajor)
|
|
|
|
refM2(j,i) = Scalar(1);
|
|
|
|
else
|
|
|
|
refM2(i,j) = Scalar(1);
|
2009-01-22 02:46:04 +08:00
|
|
|
}
|
|
|
|
}
|
2009-05-04 22:25:12 +08:00
|
|
|
}
|
|
|
|
m2.finalize();
|
2009-01-22 02:46:04 +08:00
|
|
|
VERIFY(countFalseNonZero+countTrueNonZero == m2.nonZeros());
|
|
|
|
VERIFY_IS_APPROX(m2, refM2);
|
2010-11-02 21:33:33 +08:00
|
|
|
m2.prune(Scalar(1));
|
2009-01-22 02:46:04 +08:00
|
|
|
VERIFY(countTrueNonZero==m2.nonZeros());
|
|
|
|
VERIFY_IS_APPROX(m2, refM2);
|
|
|
|
}
|
2011-12-04 21:39:24 +08:00
|
|
|
|
2012-01-28 18:13:59 +08:00
|
|
|
// test setFromTriplets
|
|
|
|
{
|
|
|
|
typedef Triplet<Scalar,Index> TripletType;
|
|
|
|
std::vector<TripletType> triplets;
|
|
|
|
int ntriplets = rows*cols;
|
|
|
|
triplets.reserve(ntriplets);
|
|
|
|
DenseMatrix refMat(rows,cols);
|
|
|
|
refMat.setZero();
|
|
|
|
for(int i=0;i<ntriplets;++i)
|
|
|
|
{
|
2012-01-31 16:14:01 +08:00
|
|
|
int r = internal::random<int>(0,rows-1);
|
|
|
|
int c = internal::random<int>(0,cols-1);
|
2012-01-28 18:13:59 +08:00
|
|
|
Scalar v = internal::random<Scalar>();
|
2012-01-31 16:14:01 +08:00
|
|
|
triplets.push_back(TripletType(r,c,v));
|
|
|
|
refMat(r,c) += v;
|
2012-01-28 18:13:59 +08:00
|
|
|
}
|
|
|
|
SparseMatrixType m(rows,cols);
|
|
|
|
m.setFromTriplets(triplets.begin(), triplets.end());
|
|
|
|
VERIFY_IS_APPROX(m, refMat);
|
|
|
|
}
|
|
|
|
|
2011-12-04 21:39:24 +08:00
|
|
|
// test triangularView
|
|
|
|
{
|
|
|
|
DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
|
|
|
|
SparseMatrixType m2(rows, rows), m3(rows, rows);
|
|
|
|
initSparse<Scalar>(density, refMat2, m2);
|
|
|
|
refMat3 = refMat2.template triangularView<Lower>();
|
|
|
|
m3 = m2.template triangularView<Lower>();
|
|
|
|
VERIFY_IS_APPROX(m3, refMat3);
|
|
|
|
|
|
|
|
refMat3 = refMat2.template triangularView<Upper>();
|
|
|
|
m3 = m2.template triangularView<Upper>();
|
|
|
|
VERIFY_IS_APPROX(m3, refMat3);
|
|
|
|
|
|
|
|
refMat3 = refMat2.template triangularView<UnitUpper>();
|
|
|
|
m3 = m2.template triangularView<UnitUpper>();
|
|
|
|
VERIFY_IS_APPROX(m3, refMat3);
|
|
|
|
|
|
|
|
refMat3 = refMat2.template triangularView<UnitLower>();
|
|
|
|
m3 = m2.template triangularView<UnitLower>();
|
|
|
|
VERIFY_IS_APPROX(m3, refMat3);
|
|
|
|
}
|
2010-07-22 21:57:01 +08:00
|
|
|
|
2010-11-15 21:14:05 +08:00
|
|
|
// test selfadjointView
|
2011-12-04 21:39:24 +08:00
|
|
|
if(!SparseMatrixType::IsRowMajor)
|
2010-11-15 21:14:05 +08:00
|
|
|
{
|
|
|
|
DenseMatrix refMat2(rows, rows), refMat3(rows, rows);
|
|
|
|
SparseMatrixType m2(rows, rows), m3(rows, rows);
|
|
|
|
initSparse<Scalar>(density, refMat2, m2);
|
|
|
|
refMat3 = refMat2.template selfadjointView<Lower>();
|
|
|
|
m3 = m2.template selfadjointView<Lower>();
|
|
|
|
VERIFY_IS_APPROX(m3, refMat3);
|
|
|
|
}
|
|
|
|
|
2010-07-22 21:57:01 +08:00
|
|
|
// test sparseView
|
|
|
|
{
|
|
|
|
DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|
|
|
SparseMatrixType m2(rows, rows);
|
2010-11-15 21:14:05 +08:00
|
|
|
initSparse<Scalar>(density, refMat2, m2);
|
2010-07-22 21:57:01 +08:00
|
|
|
VERIFY_IS_APPROX(m2.eval(), refMat2.sparseView().eval());
|
|
|
|
}
|
2011-12-05 04:49:21 +08:00
|
|
|
|
|
|
|
// test diagonal
|
|
|
|
{
|
|
|
|
DenseMatrix refMat2 = DenseMatrix::Zero(rows, rows);
|
|
|
|
SparseMatrixType m2(rows, rows);
|
|
|
|
initSparse<Scalar>(density, refMat2, m2);
|
|
|
|
VERIFY_IS_APPROX(m2.diagonal(), refMat2.diagonal().eval());
|
|
|
|
}
|
2012-07-19 06:07:06 +08:00
|
|
|
|
|
|
|
// test conservative resize
|
|
|
|
{
|
|
|
|
std::vector< std::pair<int,int> > inc;
|
|
|
|
inc.push_back(std::pair<int,int>(-3,-2));
|
|
|
|
inc.push_back(std::pair<int,int>(0,0));
|
|
|
|
inc.push_back(std::pair<int,int>(3,2));
|
|
|
|
inc.push_back(std::pair<int,int>(3,0));
|
|
|
|
inc.push_back(std::pair<int,int>(0,3));
|
|
|
|
|
|
|
|
for(size_t i = 0; i< inc.size(); i++) {
|
|
|
|
int incRows = inc[i].first;
|
|
|
|
int incCols = inc[i].second;
|
|
|
|
SparseMatrixType m1(rows, cols);
|
|
|
|
DenseMatrix refMat1 = DenseMatrix::Zero(rows, cols);
|
|
|
|
initSparse<Scalar>(density, refMat1, m1);
|
|
|
|
|
|
|
|
m1.conservativeResize(rows+incRows, cols+incCols);
|
|
|
|
refMat1.conservativeResize(rows+incRows, cols+incCols);
|
|
|
|
if (incRows > 0) refMat1.bottomRows(incRows).setZero();
|
|
|
|
if (incCols > 0) refMat1.rightCols(incCols).setZero();
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(m1, refMat1);
|
|
|
|
|
|
|
|
// Insert new values
|
|
|
|
if (incRows > 0)
|
|
|
|
m1.insert(refMat1.rows()-1, 0) = refMat1(refMat1.rows()-1, 0) = 1;
|
|
|
|
if (incCols > 0)
|
|
|
|
m1.insert(0, refMat1.cols()-1) = refMat1(0, refMat1.cols()-1) = 1;
|
|
|
|
|
|
|
|
VERIFY_IS_APPROX(m1, refMat1);
|
|
|
|
|
|
|
|
|
|
|
|
}
|
|
|
|
}
|
2008-08-22 02:40:56 +08:00
|
|
|
}
|
|
|
|
|
2008-11-05 21:47:55 +08:00
|
|
|
void test_sparse_basic()
|
2008-08-22 02:40:56 +08:00
|
|
|
{
|
2008-10-21 01:03:09 +08:00
|
|
|
for(int i = 0; i < g_repeat; i++) {
|
2011-06-07 17:28:16 +08:00
|
|
|
int s = Eigen::internal::random<int>(1,50);
|
|
|
|
CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(8, 8)) ));
|
2011-12-04 21:39:24 +08:00
|
|
|
CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, ColMajor>(s, s)) ));
|
|
|
|
CALL_SUBTEST_2(( sparse_basic(SparseMatrix<std::complex<double>, RowMajor>(s, s)) ));
|
2011-06-07 17:28:16 +08:00
|
|
|
CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double>(s, s)) ));
|
|
|
|
CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,ColMajor,long int>(s, s)) ));
|
2011-12-04 21:39:24 +08:00
|
|
|
CALL_SUBTEST_1(( sparse_basic(SparseMatrix<double,RowMajor,long int>(s, s)) ));
|
2008-10-21 01:03:09 +08:00
|
|
|
}
|
2008-08-22 01:02:47 +08:00
|
|
|
}
|