mirror of
https://gitlab.com/libeigen/eigen.git
synced 2024-12-15 07:10:37 +08:00
Add SparseVector::conservativeResize() method.
This commit is contained in:
parent
e3a15a03a4
commit
869b4443ac
@ -222,6 +222,24 @@ class SparseVector
|
||||
m_data.clear();
|
||||
}
|
||||
|
||||
/** Resizes the sparse vector to \a newSize, while leaving old values untouched.
|
||||
*
|
||||
* If the size of the vector is decreased, then the storage of the out-of bounds coefficients is kept and reserved.
|
||||
* Call .data().squeeze() to free extra memory.
|
||||
*
|
||||
* \sa reserve(), setZero()
|
||||
*/
|
||||
void conservativeResize(Index newSize)
|
||||
{
|
||||
if (newSize < m_size)
|
||||
{
|
||||
Index i = 0;
|
||||
while (i<m_data.size() && m_data.index(i)<newSize) ++i;
|
||||
m_data.resize(i);
|
||||
}
|
||||
m_size = newSize;
|
||||
}
|
||||
|
||||
void resizeNonZeros(Index size) { m_data.resize(size); }
|
||||
|
||||
inline SparseVector() : m_size(0) { check_template_parameters(); resize(0); }
|
||||
|
@ -9,14 +9,14 @@
|
||||
|
||||
#include "sparse.h"
|
||||
|
||||
template<typename Scalar,typename Index> void sparse_vector(int rows, int cols)
|
||||
template<typename Scalar,typename StorageIndex> void sparse_vector(int rows, int cols)
|
||||
{
|
||||
double densityMat = (std::max)(8./(rows*cols), 0.01);
|
||||
double densityVec = (std::max)(8./float(rows), 0.1);
|
||||
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
||||
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
||||
typedef SparseVector<Scalar,0,Index> SparseVectorType;
|
||||
typedef SparseMatrix<Scalar,0,Index> SparseMatrixType;
|
||||
typedef SparseVector<Scalar,0,StorageIndex> SparseVectorType;
|
||||
typedef SparseMatrix<Scalar,0,StorageIndex> SparseMatrixType;
|
||||
Scalar eps = 1e-6;
|
||||
|
||||
SparseMatrixType m1(rows,rows);
|
||||
@ -87,8 +87,10 @@ template<typename Scalar,typename Index> void sparse_vector(int rows, int cols)
|
||||
|
||||
VERIFY_IS_APPROX(m1*v2, refM1*refV2);
|
||||
VERIFY_IS_APPROX(v1.dot(m1*v2), refV1.dot(refM1*refV2));
|
||||
int i = internal::random<int>(0,rows-1);
|
||||
VERIFY_IS_APPROX(v1.dot(m1.col(i)), refV1.dot(refM1.col(i)));
|
||||
{
|
||||
int i = internal::random<int>(0,rows-1);
|
||||
VERIFY_IS_APPROX(v1.dot(m1.col(i)), refV1.dot(refM1.col(i)));
|
||||
}
|
||||
|
||||
|
||||
VERIFY_IS_APPROX(v1.squaredNorm(), refV1.squaredNorm());
|
||||
@ -111,15 +113,51 @@ template<typename Scalar,typename Index> void sparse_vector(int rows, int cols)
|
||||
VERIFY_IS_APPROX(refV3 = v1.transpose(),v1.toDense());
|
||||
VERIFY_IS_APPROX(DenseVector(v1),v1.toDense());
|
||||
|
||||
// test conservative resize
|
||||
{
|
||||
std::vector<StorageIndex> inc;
|
||||
if(rows > 3)
|
||||
inc.push_back(-3);
|
||||
inc.push_back(0);
|
||||
inc.push_back(3);
|
||||
inc.push_back(1);
|
||||
inc.push_back(10);
|
||||
|
||||
for(std::size_t i = 0; i< inc.size(); i++) {
|
||||
StorageIndex incRows = inc[i];
|
||||
SparseVectorType vec1(rows);
|
||||
DenseVector refVec1 = DenseVector::Zero(rows);
|
||||
initSparse<Scalar>(densityVec, refVec1, vec1);
|
||||
|
||||
vec1.conservativeResize(rows+incRows);
|
||||
refVec1.conservativeResize(rows+incRows);
|
||||
if (incRows > 0) refVec1.tail(incRows).setZero();
|
||||
|
||||
VERIFY_IS_APPROX(vec1, refVec1);
|
||||
|
||||
// Insert new values
|
||||
if (incRows > 0)
|
||||
vec1.insert(vec1.rows()-1) = refVec1(refVec1.rows()-1) = 1;
|
||||
|
||||
VERIFY_IS_APPROX(vec1, refVec1);
|
||||
}
|
||||
}
|
||||
|
||||
}
|
||||
|
||||
void test_sparse_vector()
|
||||
{
|
||||
for(int i = 0; i < g_repeat; i++) {
|
||||
int r = Eigen::internal::random<int>(1,500), c = Eigen::internal::random<int>(1,500);
|
||||
if(Eigen::internal::random<int>(0,4) == 0) {
|
||||
r = c; // check square matrices in 25% of tries
|
||||
}
|
||||
EIGEN_UNUSED_VARIABLE(r+c);
|
||||
|
||||
CALL_SUBTEST_1(( sparse_vector<double,int>(8, 8) ));
|
||||
CALL_SUBTEST_2(( sparse_vector<std::complex<double>, int>(16, 16) ));
|
||||
CALL_SUBTEST_1(( sparse_vector<double,long int>(299, 535) ));
|
||||
CALL_SUBTEST_1(( sparse_vector<double,short>(299, 535) ));
|
||||
CALL_SUBTEST_2(( sparse_vector<std::complex<double>, int>(r, c) ));
|
||||
CALL_SUBTEST_1(( sparse_vector<double,long int>(r, c) ));
|
||||
CALL_SUBTEST_1(( sparse_vector<double,short>(r, c) ));
|
||||
}
|
||||
}
|
||||
|
||||
|
Loading…
Reference in New Issue
Block a user