add unit tests for sparse LU and fix a couple of warnings

This commit is contained in:
Gael Guennebaud 2008-10-20 11:37:45 +00:00
parent fa27cd1ed0
commit e1c50a3cb1
5 changed files with 61 additions and 37 deletions

View File

@ -124,7 +124,7 @@ void AmbiVector<Scalar>::nonZeros() const
template<typename Scalar>
void AmbiVector<Scalar>::init(RealScalar estimatedDensity)
{
if (m_mode = estimatedDensity>0.1)
if (estimatedDensity>0.1)
init(IsDense);
else
init(IsSparse);

View File

@ -43,7 +43,7 @@ struct ei_traits<LinkedVectorMatrix<_Scalar,_Flags> >
template<typename Element, int ChunkSize = 8>
struct LinkedVectorChunk
{
LinkedVectorChunk() : size(0), next(0), prev(0) {}
LinkedVectorChunk() : next(0), prev(0), size(0) {}
Element data[ChunkSize];
LinkedVectorChunk* next;
LinkedVectorChunk* prev;
@ -141,7 +141,7 @@ class LinkedVectorMatrix
inline void startFill(int reserveSize = 1000)
{
clear();
for (int i=0; i<m_data.size(); ++i)
for (unsigned int i=0; i<m_data.size(); ++i)
m_ends[i] = m_data[i] = 0;
}
@ -202,7 +202,7 @@ class LinkedVectorMatrix
void clear()
{
for (int i=0; i<m_data.size(); ++i)
for (unsigned int i=0; i<m_data.size(); ++i)
{
VectorChunk* el = m_data[i];
while (el)
@ -224,7 +224,7 @@ class LinkedVectorMatrix
clear();
m_data.resize(outers);
m_ends.resize(outers);
for (int i=0; i<m_data.size(); ++i)
for (unsigned int i=0; i<m_data.size(); ++i)
m_ends[i] = m_data[i] = 0;
}
m_innerSize = inners;

View File

@ -64,11 +64,8 @@ struct SluMatrix : SuperMatrix
SluMatrix(const SluMatrix& other)
: SuperMatrix(other)
{
Store = &nnz;
nnz = other.nnz;
values = other.values;
innerInd = other.innerInd;
outerInd = other.outerInd;
Store = &storage;
storage = other.storage;
}
struct
@ -77,13 +74,13 @@ struct SluMatrix : SuperMatrix
void *values;
int *innerInd;
int *outerInd;
};
} storage;
void setStorageType(Stype_t t)
{
Stype = t;
if (t==SLU_NC || t==SLU_NR || t==SLU_DN)
Store = &nnz;
Store = &storage;
else
{
ei_assert(false && "storage type not supported");
@ -131,8 +128,8 @@ struct SluMatrixMapHelper<Matrix<Scalar,Rows,Cols,StorageOrder,MRows,MCols> >
res.nrow = mat.rows();
res.ncol = mat.cols();
res.lda = mat.stride();
res.values = mat.data();
res.storage.lda = mat.stride();
res.storage.values = mat.data();
}
};
@ -157,10 +154,10 @@ struct SluMatrixMapHelper<SparseMatrix<Scalar,Flags> >
res.Mtype = SLU_GE;
res.nnz = mat.nonZeros();
res.values = mat._valuePtr();
res.innerInd = mat._innerIndexPtr();
res.outerInd = mat._outerIndexPtr();
res.storage.nnz = mat.nonZeros();
res.storage.values = mat._valuePtr();
res.storage.innerInd = mat._innerIndexPtr();
res.storage.outerInd = mat._outerIndexPtr();
res.setScalarType<Scalar>();
@ -196,11 +193,11 @@ SparseMatrix<Scalar,Flags> SparseMatrix<Scalar,Flags>::Map(SluMatrix& sluMat)
res.m_innerSize = sluMat.nrow;
res.m_outerSize = sluMat.ncol;
}
res.m_outerIndex = sluMat.outerInd;
res.m_outerIndex = sluMat.storage.outerInd;
SparseArray<Scalar> data = SparseArray<Scalar>::Map(
sluMat.innerInd,
reinterpret_cast<Scalar*>(sluMat.values),
sluMat.outerInd[res.m_outerSize]);
sluMat.storage.innerInd,
reinterpret_cast<Scalar*>(sluMat.storage.values),
sluMat.storage.outerInd[res.m_outerSize]);
res.m_data.swap(data);
res.markAsRValue();
return res;
@ -295,9 +292,9 @@ void SparseLU<MatrixType,SuperLU>::compute(const MatrixType& a)
m_sluB.setStorageType(SLU_DN);
m_sluB.setScalarType<Scalar>();
m_sluB.Mtype = SLU_GE;
m_sluB.values = 0;
m_sluB.storage.values = 0;
m_sluB.nrow = m_sluB.ncol = 0;
m_sluB.lda = size;
m_sluB.storage.lda = size;
m_sluX = m_sluB;
StatInit(&m_sluStat);

View File

@ -23,15 +23,19 @@ if(CHOLMOD_FOUND)
set(EXTERNAL_LIBS ${EXTERNAL_LIBS} ${CHOLMOD_LIBRARIES})
endif(CHOLMOD_FOUND)
# find_package(Umfpack)
# if(UMFPACK_FOUND)
# add_definitions("-DEIGEN_UMFPACK_SUPPORT")
# endif(UMFPACK_FOUND)
#
# find_package(SuperLU)
# if(SUPERLU_FOUND)
# add_definitions("-DEIGEN_SUPERLU_SUPPORT")
# endif(SUPERLU_FOUND)
find_package(Umfpack)
if(UMFPACK_FOUND)
add_definitions("-DEIGEN_UMFPACK_SUPPORT")
include_directories(${UMFPACK_INCLUDES})
set(EXTERNAL_LIBS ${EXTERNAL_LIBS} ${UMFPACK_LIBRARIES})
endif(UMFPACK_FOUND)
find_package(SuperLU)
if(SUPERLU_FOUND)
add_definitions("-DEIGEN_SUPERLU_SUPPORT")
include_directories(${SUPERLU_INCLUDES})
set(EXTERNAL_LIBS ${EXTERNAL_LIBS} ${SUPERLU_LIBRARIES})
endif(SUPERLU_FOUND)
if(CMAKE_COMPILER_IS_GNUCXX)
if(CMAKE_SYSTEM_NAME MATCHES Linux)

View File

@ -99,7 +99,7 @@ template<typename Scalar> void sparse(int rows, int cols)
refMat.coeffRef(nonzeroCoords[0].x(), nonzeroCoords[0].y()) = Scalar(5);
VERIFY_IS_APPROX(m, refMat);
#if 0
// test InnerIterators and Block expressions
for(int j=0; j<cols; j++)
{
@ -217,7 +217,6 @@ template<typename Scalar> void sparse(int rows, int cols)
// TODO test row major
}
#endif
// test LLT
{
@ -232,8 +231,6 @@ template<typename Scalar> void sparse(int rows, int cols)
refMat2.diagonal() *= 0.5;
refMat2.llt().solve(b, &refX);
// std::cerr << refMat2 << "\n\n" << refMat2.llt().matrixL() << "\n\n";
// std::cerr << m2 << "\n\n";
typedef SparseMatrix<Scalar,Lower|SelfAdjoint> SparseSelfAdjointMatrix;
x = b;
SparseLLT<SparseSelfAdjointMatrix> (m2).solveInPlace(x);
@ -256,6 +253,32 @@ template<typename Scalar> void sparse(int rows, int cols)
#endif
}
// test LU
{
SparseMatrix<Scalar> m2(rows, cols);
DenseMatrix refMat2(rows, cols);
DenseVector b = DenseVector::Random(cols);
DenseVector refX(cols), x(cols);
initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag, &zeroCoords, &nonzeroCoords);
refMat2.lu().solve(b, &refX);
// x.setZero();
// SparseLU<SparseMatrix<Scalar> > (m2).solve(b,&x);
// VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: default");
#ifdef EIGEN_SUPERLU_SUPPORT
x.setZero();
SparseLU<SparseMatrix<Scalar>,SuperLU>(m2).solve(b,&x);
VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: SuperLU");
#endif
#ifdef EIGEN_UMFPACK_SUPPORT
x.setZero();
SparseLU<SparseMatrix<Scalar>,UmfPack>(m2).solve(b,&x);
VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: umfpack");
#endif
}
}
void test_sparse()