2011-10-11 17:32:26 +08:00
|
|
|
// This file is part of Eigen, a lightweight C++ template library
|
|
|
|
// for linear algebra.
|
|
|
|
//
|
|
|
|
// Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
|
|
|
|
//
|
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/.
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
#include "sparse.h"
|
2011-11-12 21:11:27 +08:00
|
|
|
#include <Eigen/SparseCore>
|
2015-06-05 20:33:37 +08:00
|
|
|
#include <sstream>
|
2011-10-11 17:32:26 +08:00
|
|
|
|
2016-07-04 23:19:38 +08:00
|
|
|
template<typename Solver, typename Rhs, typename Guess,typename Result>
|
|
|
|
void solve_with_guess(IterativeSolverBase<Solver>& solver, const MatrixBase<Rhs>& b, const Guess& g, Result &x) {
|
2016-07-04 23:47:47 +08:00
|
|
|
if(internal::random<bool>())
|
|
|
|
{
|
|
|
|
// With a temporary through evaluator<SolveWithGuess>
|
|
|
|
x = solver.derived().solveWithGuess(b,g) + Result::Zero(x.rows(), x.cols());
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
// direct evaluation within x through Assignment<Result,SolveWithGuess>
|
|
|
|
x = solver.derived().solveWithGuess(b.derived(),g);
|
|
|
|
}
|
2016-07-04 23:19:38 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver, typename Rhs, typename Guess,typename Result>
|
|
|
|
void solve_with_guess(SparseSolverBase<Solver>& solver, const MatrixBase<Rhs>& b, const Guess& , Result& x) {
|
2016-07-04 23:47:47 +08:00
|
|
|
if(internal::random<bool>())
|
|
|
|
x = solver.derived().solve(b) + Result::Zero(x.rows(), x.cols());
|
|
|
|
else
|
|
|
|
x = solver.derived().solve(b);
|
2016-07-04 23:19:38 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver, typename Rhs, typename Guess,typename Result>
|
|
|
|
void solve_with_guess(SparseSolverBase<Solver>& solver, const SparseMatrixBase<Rhs>& b, const Guess& , Result& x) {
|
|
|
|
x = solver.derived().solve(b);
|
|
|
|
}
|
|
|
|
|
2011-10-11 17:32:26 +08:00
|
|
|
template<typename Solver, typename Rhs, typename DenseMat, typename DenseRhs>
|
|
|
|
void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const DenseMat& dA, const DenseRhs& db)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
2014-12-05 05:48:53 +08:00
|
|
|
typedef typename Mat::StorageIndex StorageIndex;
|
2011-10-11 17:32:26 +08:00
|
|
|
|
2015-03-04 16:34:27 +08:00
|
|
|
DenseRhs refX = dA.householderQr().solve(db);
|
2013-06-24 01:11:32 +08:00
|
|
|
{
|
2015-03-04 16:34:27 +08:00
|
|
|
Rhs x(A.cols(), b.cols());
|
2013-06-24 01:11:32 +08:00
|
|
|
Rhs oldb = b;
|
2011-10-11 17:32:26 +08:00
|
|
|
|
2013-06-24 01:11:32 +08:00
|
|
|
solver.compute(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2015-06-05 20:33:37 +08:00
|
|
|
std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
|
|
|
|
VERIFY(solver.info() == Success);
|
2013-06-24 01:11:32 +08:00
|
|
|
}
|
|
|
|
x = solver.solve(b);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2018-10-16 06:43:44 +08:00
|
|
|
std::cerr << "WARNING: sparse solver testing: solving failed (" << typeid(Solver).name() << ")\n";
|
|
|
|
// dump call stack:
|
|
|
|
g_test_level++;
|
|
|
|
VERIFY(solver.info() == Success);
|
|
|
|
g_test_level--;
|
2013-06-24 01:11:32 +08:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
2016-07-04 23:19:38 +08:00
|
|
|
|
|
|
|
x.setZero();
|
|
|
|
solve_with_guess(solver, b, x, x);
|
2018-10-16 06:43:44 +08:00
|
|
|
VERIFY(solver.info() == Success && "solving failed when using solve_with_guess API");
|
2016-07-04 23:19:38 +08:00
|
|
|
VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
2014-10-20 22:46:47 +08:00
|
|
|
|
2013-06-24 01:11:32 +08:00
|
|
|
x.setZero();
|
|
|
|
// test the analyze/factorize API
|
|
|
|
solver.analyzePattern(A);
|
|
|
|
solver.factorize(A);
|
2015-06-05 20:33:37 +08:00
|
|
|
VERIFY(solver.info() == Success && "factorization failed when using analyzePattern/factorize API");
|
2013-06-24 01:11:32 +08:00
|
|
|
x = solver.solve(b);
|
2015-06-05 20:33:37 +08:00
|
|
|
VERIFY(solver.info() == Success && "solving failed when using analyzePattern/factorize API");
|
2013-06-24 01:11:32 +08:00
|
|
|
VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
2014-10-20 22:46:47 +08:00
|
|
|
|
|
|
|
x.setZero();
|
|
|
|
// test with Map
|
2014-12-05 05:48:53 +08:00
|
|
|
MappedSparseMatrix<Scalar,Mat::Options,StorageIndex> Am(A.rows(), A.cols(), A.nonZeros(), const_cast<StorageIndex*>(A.outerIndexPtr()), const_cast<StorageIndex*>(A.innerIndexPtr()), const_cast<Scalar*>(A.valuePtr()));
|
2014-10-20 22:46:47 +08:00
|
|
|
solver.compute(Am);
|
2015-06-05 20:33:37 +08:00
|
|
|
VERIFY(solver.info() == Success && "factorization failed when using Map");
|
2014-10-20 22:46:47 +08:00
|
|
|
DenseRhs dx(refX);
|
|
|
|
dx.setZero();
|
|
|
|
Map<DenseRhs> xm(dx.data(), dx.rows(), dx.cols());
|
|
|
|
Map<const DenseRhs> bm(db.data(), db.rows(), db.cols());
|
|
|
|
xm = solver.solve(bm);
|
2015-06-05 20:33:37 +08:00
|
|
|
VERIFY(solver.info() == Success && "solving failed when using Map");
|
2014-10-20 22:46:47 +08:00
|
|
|
VERIFY(oldb.isApprox(bm) && "sparse solver testing: the rhs should not be modified!");
|
|
|
|
VERIFY(xm.isApprox(refX,test_precision<Scalar>()));
|
2012-02-27 20:22:38 +08:00
|
|
|
}
|
2012-06-06 15:40:01 +08:00
|
|
|
|
2015-10-26 23:16:24 +08:00
|
|
|
// if not too large, do some extra check:
|
|
|
|
if(A.rows()<2000)
|
2015-02-17 00:05:10 +08:00
|
|
|
{
|
2015-10-26 23:16:24 +08:00
|
|
|
// test initialization ctor
|
|
|
|
{
|
|
|
|
Rhs x(b.rows(), b.cols());
|
|
|
|
Solver solver2(A);
|
|
|
|
VERIFY(solver2.info() == Success);
|
|
|
|
x = solver2.solve(b);
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
|
|
|
}
|
|
|
|
|
|
|
|
// test dense Block as the result and rhs:
|
|
|
|
{
|
|
|
|
DenseRhs x(refX.rows(), refX.cols());
|
|
|
|
DenseRhs oldb(db);
|
|
|
|
x.setZero();
|
|
|
|
x.block(0,0,x.rows(),x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
|
|
|
|
VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!");
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
|
|
|
}
|
|
|
|
|
|
|
|
// test uncompressed inputs
|
|
|
|
{
|
|
|
|
Mat A2 = A;
|
|
|
|
A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().eval());
|
|
|
|
solver.compute(A2);
|
|
|
|
Rhs x = solver.solve(b);
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
|
|
|
}
|
|
|
|
|
|
|
|
// test expression as input
|
|
|
|
{
|
|
|
|
solver.compute(0.5*(A+A));
|
|
|
|
Rhs x = solver.solve(b);
|
|
|
|
VERIFY(x.isApprox(refX,test_precision<Scalar>()));
|
|
|
|
|
|
|
|
Solver solver2(0.5*(A+A));
|
|
|
|
Rhs x2 = solver2.solve(b);
|
|
|
|
VERIFY(x2.isApprox(refX,test_precision<Scalar>()));
|
|
|
|
}
|
2014-10-06 17:41:50 +08:00
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
}
|
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
template<typename Solver, typename Rhs>
|
2015-06-05 20:33:37 +08:00
|
|
|
void check_sparse_solving_real_cases(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const typename Solver::MatrixType& fullA, const Rhs& refX)
|
2012-03-29 20:32:54 +08:00
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef typename Mat::RealScalar RealScalar;
|
|
|
|
|
2015-03-04 16:34:27 +08:00
|
|
|
Rhs x(A.cols(), b.cols());
|
2015-06-05 20:33:37 +08:00
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
solver.compute(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2015-06-05 20:33:37 +08:00
|
|
|
std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
|
|
|
|
VERIFY(solver.info() == Success);
|
2012-03-29 20:32:54 +08:00
|
|
|
}
|
|
|
|
x = solver.solve(b);
|
2015-06-05 20:33:37 +08:00
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2015-06-05 20:33:37 +08:00
|
|
|
std::cerr << "WARNING | sparse solver testing, solving failed (" << typeid(Solver).name() << ")\n";
|
2012-03-29 20:32:54 +08:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
2015-06-05 20:33:37 +08:00
|
|
|
RealScalar res_error = (fullA*x-b).norm()/b.norm();
|
|
|
|
VERIFY( (res_error <= test_precision<Scalar>() ) && "sparse solver failed without noticing it");
|
|
|
|
|
|
|
|
|
|
|
|
if(refX.size() != 0 && (refX - x).norm()/refX.norm() > test_precision<Scalar>())
|
|
|
|
{
|
|
|
|
std::cerr << "WARNING | found solution is different from the provided reference one\n";
|
2012-03-29 20:32:54 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
template<typename Solver, typename DenseMat>
|
|
|
|
void check_sparse_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
|
|
|
|
solver.compute(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2015-06-05 20:33:37 +08:00
|
|
|
std::cerr << "WARNING | sparse solver testing: factorization failed (check_sparse_determinant)\n";
|
2011-10-11 17:32:26 +08:00
|
|
|
return;
|
|
|
|
}
|
|
|
|
|
|
|
|
Scalar refDet = dA.determinant();
|
|
|
|
VERIFY_IS_APPROX(refDet,solver.determinant());
|
|
|
|
}
|
2014-10-17 22:52:56 +08:00
|
|
|
template<typename Solver, typename DenseMat>
|
|
|
|
void check_sparse_abs_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA)
|
|
|
|
{
|
|
|
|
using std::abs;
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
|
|
|
|
solver.compute(A);
|
|
|
|
if (solver.info() != Success)
|
|
|
|
{
|
2015-06-05 20:33:37 +08:00
|
|
|
std::cerr << "WARNING | sparse solver testing: factorization failed (check_sparse_abs_determinant)\n";
|
2014-10-17 22:52:56 +08:00
|
|
|
return;
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
|
2014-10-17 22:52:56 +08:00
|
|
|
Scalar refDet = abs(dA.determinant());
|
|
|
|
VERIFY_IS_APPROX(refDet,solver.absDeterminant());
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
template<typename Solver, typename DenseMat>
|
|
|
|
int generate_sparse_spd_problem(Solver& , typename Solver::MatrixType& A, typename Solver::MatrixType& halfA, DenseMat& dA, int maxSize = 300)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
|
|
|
|
int size = internal::random<int>(1,maxSize);
|
|
|
|
double density = (std::max)(8./(size*size), 0.01);
|
|
|
|
|
|
|
|
Mat M(size, size);
|
|
|
|
DenseMatrix dM(size, size);
|
|
|
|
|
|
|
|
initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
|
|
|
|
|
|
|
|
A = M * M.adjoint();
|
|
|
|
dA = dM * dM.adjoint();
|
|
|
|
|
|
|
|
halfA.resize(size,size);
|
2015-02-11 01:57:41 +08:00
|
|
|
if(Solver::UpLo==(Lower|Upper))
|
|
|
|
halfA = A;
|
|
|
|
else
|
|
|
|
halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
return size;
|
|
|
|
}
|
|
|
|
|
2012-03-29 20:53:42 +08:00
|
|
|
|
|
|
|
#ifdef TEST_REAL_CASES
|
|
|
|
template<typename Scalar>
|
|
|
|
inline std::string get_matrixfolder()
|
|
|
|
{
|
|
|
|
std::string mat_folder = TEST_REAL_CASES;
|
|
|
|
if( internal::is_same<Scalar, std::complex<float> >::value || internal::is_same<Scalar, std::complex<double> >::value )
|
2012-08-03 19:05:45 +08:00
|
|
|
mat_folder = mat_folder + static_cast<std::string>("/complex/");
|
2012-03-29 20:53:42 +08:00
|
|
|
else
|
2012-08-03 19:05:45 +08:00
|
|
|
mat_folder = mat_folder + static_cast<std::string>("/real/");
|
2012-03-29 20:53:42 +08:00
|
|
|
return mat_folder;
|
|
|
|
}
|
2015-06-05 20:33:37 +08:00
|
|
|
std::string sym_to_string(int sym)
|
|
|
|
{
|
|
|
|
if(sym==Symmetric) return "Symmetric ";
|
|
|
|
if(sym==SPD) return "SPD ";
|
|
|
|
return "";
|
|
|
|
}
|
|
|
|
template<typename Derived>
|
|
|
|
std::string solver_stats(const IterativeSolverBase<Derived> &solver)
|
|
|
|
{
|
|
|
|
std::stringstream ss;
|
|
|
|
ss << solver.iterations() << " iters, error: " << solver.error();
|
|
|
|
return ss.str();
|
|
|
|
}
|
|
|
|
template<typename Derived>
|
|
|
|
std::string solver_stats(const SparseSolverBase<Derived> &/*solver*/)
|
|
|
|
{
|
|
|
|
return "";
|
|
|
|
}
|
2012-03-29 20:53:42 +08:00
|
|
|
#endif
|
|
|
|
|
2018-04-14 02:36:58 +08:00
|
|
|
template<typename Solver> void check_sparse_spd_solving(Solver& solver, int maxSize = (std::min)(300,EIGEN_TEST_MAX_SIZE), int maxRealWorldSize = 100000)
|
2011-10-11 17:32:26 +08:00
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
2015-06-05 20:33:37 +08:00
|
|
|
typedef typename Mat::StorageIndex StorageIndex;
|
|
|
|
typedef SparseMatrix<Scalar,ColMajor, StorageIndex> SpMat;
|
2016-11-07 03:29:57 +08:00
|
|
|
typedef SparseVector<Scalar, 0, StorageIndex> SpVec;
|
2011-10-11 17:32:26 +08:00
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
|
|
|
|
|
|
|
// generate the problem
|
|
|
|
Mat A, halfA;
|
|
|
|
DenseMatrix dA;
|
2012-03-29 20:32:54 +08:00
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
2015-06-09 15:29:53 +08:00
|
|
|
int size = generate_sparse_spd_problem(solver, A, halfA, dA, maxSize);
|
2013-03-20 23:15:18 +08:00
|
|
|
|
|
|
|
// generate the right hand sides
|
|
|
|
int rhsCols = internal::random<int>(1,16);
|
|
|
|
double density = (std::max)(8./(size*rhsCols), 0.1);
|
|
|
|
SpMat B(size,rhsCols);
|
|
|
|
DenseVector b = DenseVector::Random(size);
|
|
|
|
DenseMatrix dB(size,rhsCols);
|
|
|
|
initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
|
2016-11-07 03:29:57 +08:00
|
|
|
SpVec c = B.col(0);
|
|
|
|
DenseVector dc = dB.col(0);
|
2013-03-20 23:15:18 +08:00
|
|
|
|
2015-06-05 20:33:37 +08:00
|
|
|
CALL_SUBTEST( check_sparse_solving(solver, A, b, dA, b) );
|
|
|
|
CALL_SUBTEST( check_sparse_solving(solver, halfA, b, dA, b) );
|
|
|
|
CALL_SUBTEST( check_sparse_solving(solver, A, dB, dA, dB) );
|
|
|
|
CALL_SUBTEST( check_sparse_solving(solver, halfA, dB, dA, dB) );
|
|
|
|
CALL_SUBTEST( check_sparse_solving(solver, A, B, dA, dB) );
|
|
|
|
CALL_SUBTEST( check_sparse_solving(solver, halfA, B, dA, dB) );
|
2016-11-07 03:29:57 +08:00
|
|
|
CALL_SUBTEST( check_sparse_solving(solver, A, c, dA, dc) );
|
|
|
|
CALL_SUBTEST( check_sparse_solving(solver, halfA, c, dA, dc) );
|
2013-03-20 23:15:18 +08:00
|
|
|
|
|
|
|
// check only once
|
|
|
|
if(i==0)
|
|
|
|
{
|
|
|
|
b = DenseVector::Zero(size);
|
|
|
|
check_sparse_solving(solver, A, b, dA, b);
|
|
|
|
}
|
2012-03-29 20:32:54 +08:00
|
|
|
}
|
2013-03-20 23:15:18 +08:00
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
// First, get the folder
|
2015-06-05 20:33:37 +08:00
|
|
|
#ifdef TEST_REAL_CASES
|
|
|
|
// Test real problems with double precision only
|
|
|
|
if (internal::is_same<typename NumTraits<Scalar>::Real, double>::value)
|
2012-03-29 20:32:54 +08:00
|
|
|
{
|
2015-06-05 20:33:37 +08:00
|
|
|
std::string mat_folder = get_matrixfolder<Scalar>();
|
|
|
|
MatrixMarketIterator<Scalar> it(mat_folder);
|
|
|
|
for (; it; ++it)
|
|
|
|
{
|
|
|
|
if (it.sym() == SPD){
|
|
|
|
A = it.matrix();
|
2015-06-09 15:29:53 +08:00
|
|
|
if(A.diagonal().size() <= maxRealWorldSize)
|
|
|
|
{
|
|
|
|
DenseVector b = it.rhs();
|
|
|
|
DenseVector refX = it.refX();
|
|
|
|
PermutationMatrix<Dynamic, Dynamic, StorageIndex> pnull;
|
2015-06-25 23:08:58 +08:00
|
|
|
halfA.resize(A.rows(), A.cols());
|
2015-06-09 15:29:53 +08:00
|
|
|
if(Solver::UpLo == (Lower|Upper))
|
|
|
|
halfA = A;
|
|
|
|
else
|
|
|
|
halfA.template selfadjointView<Solver::UpLo>() = A.template triangularView<Eigen::Lower>().twistedBy(pnull);
|
|
|
|
|
|
|
|
std::cout << "INFO | Testing " << sym_to_string(it.sym()) << "sparse problem " << it.matname()
|
|
|
|
<< " (" << A.rows() << "x" << A.cols() << ") using " << typeid(Solver).name() << "..." << std::endl;
|
|
|
|
CALL_SUBTEST( check_sparse_solving_real_cases(solver, A, b, A, refX) );
|
|
|
|
std::string stats = solver_stats(solver);
|
|
|
|
if(stats.size()>0)
|
|
|
|
std::cout << "INFO | " << stats << std::endl;
|
|
|
|
CALL_SUBTEST( check_sparse_solving_real_cases(solver, halfA, b, A, refX) );
|
|
|
|
}
|
2015-06-05 20:33:37 +08:00
|
|
|
else
|
2015-06-09 15:29:53 +08:00
|
|
|
{
|
|
|
|
std::cout << "INFO | Skip sparse problem \"" << it.matname() << "\" (too large)" << std::endl;
|
|
|
|
}
|
2015-06-05 20:33:37 +08:00
|
|
|
}
|
2012-03-29 20:32:54 +08:00
|
|
|
}
|
|
|
|
}
|
2015-06-09 21:17:21 +08:00
|
|
|
#else
|
|
|
|
EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
|
2012-03-29 20:32:54 +08:00
|
|
|
#endif
|
2011-10-11 17:32:26 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_spd_determinant(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
|
|
|
|
// generate the problem
|
|
|
|
Mat A, halfA;
|
|
|
|
DenseMatrix dA;
|
|
|
|
generate_sparse_spd_problem(solver, A, halfA, dA, 30);
|
2012-03-29 20:32:54 +08:00
|
|
|
|
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
|
|
|
check_sparse_determinant(solver, A, dA);
|
|
|
|
check_sparse_determinant(solver, halfA, dA );
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver, typename DenseMat>
|
2015-03-09 20:54:05 +08:00
|
|
|
Index generate_sparse_square_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300, int options = ForceNonZeroDiag)
|
2011-10-11 17:32:26 +08:00
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
|
2015-03-09 20:54:05 +08:00
|
|
|
Index size = internal::random<int>(1,maxSize);
|
2011-10-11 17:32:26 +08:00
|
|
|
double density = (std::max)(8./(size*size), 0.01);
|
|
|
|
|
|
|
|
A.resize(size,size);
|
|
|
|
dA.resize(size,size);
|
|
|
|
|
2015-02-17 02:09:48 +08:00
|
|
|
initSparse<Scalar>(density, dA, A, options);
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
return size;
|
|
|
|
}
|
|
|
|
|
2015-07-27 02:30:30 +08:00
|
|
|
|
|
|
|
struct prune_column {
|
|
|
|
Index m_col;
|
|
|
|
prune_column(Index col) : m_col(col) {}
|
|
|
|
template<class Scalar>
|
|
|
|
bool operator()(Index, Index col, const Scalar&) const {
|
|
|
|
return col != m_col;
|
|
|
|
}
|
|
|
|
};
|
|
|
|
|
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_square_solving(Solver& solver, int maxSize = 300, int maxRealWorldSize = 100000, bool checkDeficient = false)
|
2011-10-11 17:32:26 +08:00
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
2015-03-09 20:54:05 +08:00
|
|
|
typedef SparseMatrix<Scalar,ColMajor, typename Mat::StorageIndex> SpMat;
|
2016-11-07 03:29:57 +08:00
|
|
|
typedef SparseVector<Scalar, 0, typename Mat::StorageIndex> SpVec;
|
2011-10-11 17:32:26 +08:00
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
|
|
|
|
|
|
|
int rhsCols = internal::random<int>(1,16);
|
|
|
|
|
|
|
|
Mat A;
|
|
|
|
DenseMatrix dA;
|
2012-03-29 20:32:54 +08:00
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
2015-06-09 15:29:53 +08:00
|
|
|
Index size = generate_sparse_square_problem(solver, A, dA, maxSize);
|
2013-03-20 23:15:18 +08:00
|
|
|
|
|
|
|
A.makeCompressed();
|
|
|
|
DenseVector b = DenseVector::Random(size);
|
|
|
|
DenseMatrix dB(size,rhsCols);
|
|
|
|
SpMat B(size,rhsCols);
|
|
|
|
double density = (std::max)(8./(size*rhsCols), 0.1);
|
|
|
|
initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
|
|
|
|
B.makeCompressed();
|
2016-11-07 03:29:57 +08:00
|
|
|
SpVec c = B.col(0);
|
|
|
|
DenseVector dc = dB.col(0);
|
2015-06-05 20:33:37 +08:00
|
|
|
CALL_SUBTEST(check_sparse_solving(solver, A, b, dA, b));
|
|
|
|
CALL_SUBTEST(check_sparse_solving(solver, A, dB, dA, dB));
|
|
|
|
CALL_SUBTEST(check_sparse_solving(solver, A, B, dA, dB));
|
2016-11-07 03:29:57 +08:00
|
|
|
CALL_SUBTEST(check_sparse_solving(solver, A, c, dA, dc));
|
2013-03-20 23:15:18 +08:00
|
|
|
|
|
|
|
// check only once
|
|
|
|
if(i==0)
|
|
|
|
{
|
2018-10-16 06:43:44 +08:00
|
|
|
CALL_SUBTEST(b = DenseVector::Zero(size); check_sparse_solving(solver, A, b, dA, b));
|
2013-03-20 23:15:18 +08:00
|
|
|
}
|
2015-07-27 02:30:30 +08:00
|
|
|
// regression test for Bug 792 (structurally rank deficient matrices):
|
|
|
|
if(checkDeficient && size>1) {
|
2015-08-04 22:12:44 +08:00
|
|
|
Index col = internal::random<int>(0,int(size-1));
|
2015-07-27 02:30:30 +08:00
|
|
|
A.prune(prune_column(col));
|
|
|
|
solver.compute(A);
|
|
|
|
VERIFY_IS_EQUAL(solver.info(), NumericalIssue);
|
|
|
|
}
|
2012-03-29 20:32:54 +08:00
|
|
|
}
|
2013-03-20 23:15:18 +08:00
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
// First, get the folder
|
2012-03-29 20:53:42 +08:00
|
|
|
#ifdef TEST_REAL_CASES
|
2015-06-05 20:33:37 +08:00
|
|
|
// Test real problems with double precision only
|
|
|
|
if (internal::is_same<typename NumTraits<Scalar>::Real, double>::value)
|
2012-03-29 20:32:54 +08:00
|
|
|
{
|
2015-06-05 20:33:37 +08:00
|
|
|
std::string mat_folder = get_matrixfolder<Scalar>();
|
|
|
|
MatrixMarketIterator<Scalar> it(mat_folder);
|
|
|
|
for (; it; ++it)
|
|
|
|
{
|
|
|
|
A = it.matrix();
|
2015-06-09 15:29:53 +08:00
|
|
|
if(A.diagonal().size() <= maxRealWorldSize)
|
|
|
|
{
|
|
|
|
DenseVector b = it.rhs();
|
|
|
|
DenseVector refX = it.refX();
|
|
|
|
std::cout << "INFO | Testing " << sym_to_string(it.sym()) << "sparse problem " << it.matname()
|
|
|
|
<< " (" << A.rows() << "x" << A.cols() << ") using " << typeid(Solver).name() << "..." << std::endl;
|
|
|
|
CALL_SUBTEST(check_sparse_solving_real_cases(solver, A, b, A, refX));
|
|
|
|
std::string stats = solver_stats(solver);
|
|
|
|
if(stats.size()>0)
|
|
|
|
std::cout << "INFO | " << stats << std::endl;
|
|
|
|
}
|
|
|
|
else
|
|
|
|
{
|
|
|
|
std::cout << "INFO | SKIP sparse problem \"" << it.matname() << "\" (too large)" << std::endl;
|
|
|
|
}
|
2015-06-05 20:33:37 +08:00
|
|
|
}
|
2012-03-29 20:32:54 +08:00
|
|
|
}
|
2015-06-09 21:17:21 +08:00
|
|
|
#else
|
|
|
|
EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
|
2012-03-29 20:32:54 +08:00
|
|
|
#endif
|
2011-10-11 17:32:26 +08:00
|
|
|
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_square_determinant(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
2015-02-17 02:09:48 +08:00
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
2015-02-17 02:09:48 +08:00
|
|
|
// generate the problem
|
|
|
|
Mat A;
|
|
|
|
DenseMatrix dA;
|
|
|
|
|
|
|
|
int size = internal::random<int>(1,30);
|
|
|
|
dA.setRandom(size,size);
|
|
|
|
|
|
|
|
dA = (dA.array().abs()<0.3).select(0,dA);
|
|
|
|
dA.diagonal() = (dA.diagonal().array()==0).select(1,dA.diagonal());
|
|
|
|
A = dA.sparseView();
|
|
|
|
A.makeCompressed();
|
|
|
|
|
2012-03-29 20:32:54 +08:00
|
|
|
check_sparse_determinant(solver, A, dA);
|
|
|
|
}
|
2011-10-11 17:32:26 +08:00
|
|
|
}
|
2014-10-17 22:52:56 +08:00
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_square_abs_determinant(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
|
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
2015-02-17 02:09:48 +08:00
|
|
|
// generate the problem
|
|
|
|
Mat A;
|
|
|
|
DenseMatrix dA;
|
|
|
|
generate_sparse_square_problem(solver, A, dA, 30);
|
|
|
|
A.makeCompressed();
|
2014-10-17 22:52:56 +08:00
|
|
|
check_sparse_abs_determinant(solver, A, dA);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
|
2015-03-04 16:34:27 +08:00
|
|
|
template<typename Solver, typename DenseMat>
|
|
|
|
void generate_sparse_leastsquare_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300, int options = ForceNonZeroDiag)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
|
|
|
|
|
|
|
int rows = internal::random<int>(1,maxSize);
|
|
|
|
int cols = internal::random<int>(1,rows);
|
|
|
|
double density = (std::max)(8./(rows*cols), 0.01);
|
|
|
|
|
|
|
|
A.resize(rows,cols);
|
|
|
|
dA.resize(rows,cols);
|
|
|
|
|
|
|
|
initSparse<Scalar>(density, dA, A, options);
|
|
|
|
}
|
|
|
|
|
|
|
|
template<typename Solver> void check_sparse_leastsquare_solving(Solver& solver)
|
|
|
|
{
|
|
|
|
typedef typename Solver::MatrixType Mat;
|
|
|
|
typedef typename Mat::Scalar Scalar;
|
2015-03-09 20:54:05 +08:00
|
|
|
typedef SparseMatrix<Scalar,ColMajor, typename Mat::StorageIndex> SpMat;
|
2015-03-04 16:34:27 +08:00
|
|
|
typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
|
|
|
|
typedef Matrix<Scalar,Dynamic,1> DenseVector;
|
|
|
|
|
|
|
|
int rhsCols = internal::random<int>(1,16);
|
|
|
|
|
|
|
|
Mat A;
|
|
|
|
DenseMatrix dA;
|
|
|
|
for (int i = 0; i < g_repeat; i++) {
|
|
|
|
generate_sparse_leastsquare_problem(solver, A, dA);
|
|
|
|
|
|
|
|
A.makeCompressed();
|
|
|
|
DenseVector b = DenseVector::Random(A.rows());
|
|
|
|
DenseMatrix dB(A.rows(),rhsCols);
|
|
|
|
SpMat B(A.rows(),rhsCols);
|
|
|
|
double density = (std::max)(8./(A.rows()*rhsCols), 0.1);
|
|
|
|
initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
|
|
|
|
B.makeCompressed();
|
|
|
|
check_sparse_solving(solver, A, b, dA, b);
|
|
|
|
check_sparse_solving(solver, A, dB, dA, dB);
|
|
|
|
check_sparse_solving(solver, A, B, dA, dB);
|
|
|
|
|
|
|
|
// check only once
|
|
|
|
if(i==0)
|
|
|
|
{
|
|
|
|
b = DenseVector::Zero(A.rows());
|
|
|
|
check_sparse_solving(solver, A, b, dA, b);
|
|
|
|
}
|
|
|
|
}
|
|
|
|
}
|