11 #include <Eigen/SparseCore> 14 template<
typename Solver,
typename Rhs,
typename Guess,
typename Result>
16 if(internal::random<bool>())
19 x = solver.
derived().solveWithGuess(b,g) + Result::Zero(x.rows(), x.cols());
24 x = solver.
derived().solveWithGuess(b.derived(),
g);
28 template<
typename Solver,
typename Rhs,
typename Guess,
typename Result>
30 if(internal::random<bool>())
31 x = solver.
derived().solve(b) + Result::Zero(x.rows(), x.cols());
36 template<
typename Solver,
typename Rhs,
typename Guess,
typename Result>
41 template<
typename Solver,
typename Rhs,
typename DenseMat,
typename DenseRhs>
46 typedef typename Mat::StorageIndex StorageIndex;
48 DenseRhs refX = dA.householderQr().solve(db);
50 Rhs x(A.cols(), b.cols());
56 std::cerr <<
"ERROR | sparse solver testing, factorization failed (" <<
typeid(Solver).
name() <<
")\n";
62 std::cerr <<
"WARNING | sparse solver testing: solving failed (" <<
typeid(Solver).
name() <<
")\n";
65 VERIFY(oldb.isApprox(b) &&
"sparse solver testing: the rhs should not be modified!");
66 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
70 VERIFY(solver.info() ==
Success &&
"solving failed when using analyzePattern/factorize API");
71 VERIFY(oldb.isApprox(b) &&
"sparse solver testing: the rhs should not be modified!");
72 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
76 solver.analyzePattern(A);
78 VERIFY(solver.info() ==
Success &&
"factorization failed when using analyzePattern/factorize API");
80 VERIFY(solver.info() ==
Success &&
"solving failed when using analyzePattern/factorize API");
81 VERIFY(oldb.isApprox(b) &&
"sparse solver testing: the rhs should not be modified!");
82 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
88 VERIFY(solver.info() ==
Success &&
"factorization failed when using Map");
93 xm = solver.solve(bm);
94 VERIFY(solver.info() ==
Success &&
"solving failed when using Map");
95 VERIFY(oldb.isApprox(bm) &&
"sparse solver testing: the rhs should not be modified!");
96 VERIFY(xm.isApprox(refX,test_precision<Scalar>()));
104 Rhs x(b.rows(), b.cols());
107 x = solver2.solve(b);
108 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
113 DenseRhs
x(refX.rows(), refX.cols());
116 x.block(0,0,
x.rows(),
x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
117 VERIFY(oldb.isApprox(db) &&
"sparse solver testing: the rhs should not be modified!");
118 VERIFY(
x.isApprox(refX,test_precision<Scalar>()));
124 A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().
eval());
126 Rhs x = solver.solve(b);
127 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
132 solver.compute(0.5*(A+A));
133 Rhs x = solver.solve(b);
134 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
136 Solver solver2(0.5*(A+A));
137 Rhs x2 = solver2.solve(b);
138 VERIFY(x2.isApprox(refX,test_precision<Scalar>()));
143 template<
typename Solver,
typename Rhs>
150 Rhs x(A.cols(), b.cols());
155 std::cerr <<
"ERROR | sparse solver testing, factorization failed (" <<
typeid(Solver).
name() <<
")\n";
162 std::cerr <<
"WARNING | sparse solver testing, solving failed (" <<
typeid(Solver).
name() <<
")\n";
166 RealScalar res_error = (fullA*x-
b).norm()/b.norm();
167 VERIFY( (res_error <= test_precision<Scalar>() ) &&
"sparse solver failed without noticing it");
170 if(refX.size() != 0 && (refX -
x).norm()/refX.norm() > test_precision<Scalar>())
172 std::cerr <<
"WARNING | found solution is different from the provided reference one\n";
176 template<
typename Solver,
typename DenseMat>
185 std::cerr <<
"WARNING | sparse solver testing: factorization failed (check_sparse_determinant)\n";
189 Scalar refDet = dA.determinant();
192 template<
typename Solver,
typename DenseMat>
202 std::cerr <<
"WARNING | sparse solver testing: factorization failed (check_sparse_abs_determinant)\n";
206 Scalar refDet =
abs(dA.determinant());
210 template<
typename Solver,
typename DenseMat>
217 int size = internal::random<int>(1,maxSize);
218 double density = (
std::max)(8./(size*size), 0.01);
221 DenseMatrix dM(size, size);
226 dA = dM * dM.adjoint();
228 halfA.resize(size,size);
232 halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
238 #ifdef TEST_REAL_CASES 239 template<
typename Scalar>
240 inline std::string get_matrixfolder()
242 std::string mat_folder = TEST_REAL_CASES;
243 if( internal::is_same<
Scalar, std::complex<float> >::
value || internal::is_same<
Scalar, std::complex<double> >::
value )
244 mat_folder = mat_folder +
static_cast<std::string
>(
"/complex/");
246 mat_folder = mat_folder +
static_cast<std::string
>(
"/real/");
249 std::string sym_to_string(
int sym)
252 if(sym==
SPD)
return "SPD ";
255 template<
typename Derived>
258 std::stringstream
ss;
262 template<
typename Derived>
273 typedef typename Mat::StorageIndex StorageIndex;
286 int rhsCols = internal::random<int>(1,16);
287 double density = (
std::max)(8./(size*rhsCols), 0.1);
288 SpMat
B(size,rhsCols);
289 DenseVector
b = DenseVector::Random(size);
290 DenseMatrix dB(size,rhsCols);
293 DenseVector dc = dB.col(0);
307 b = DenseVector::Zero(size);
313 #ifdef TEST_REAL_CASES 317 std::string mat_folder = get_matrixfolder<Scalar>();
321 if (it.sym() ==
SPD){
323 if(A.diagonal().size() <= maxRealWorldSize)
325 DenseVector
b = it.rhs();
326 DenseVector refX = it.refX();
332 halfA.template selfadjointView<Solver::UpLo>() = A.template triangularView<Eigen::Lower>().twistedBy(pnull);
334 std::cout <<
"INFO | Testing " << sym_to_string(it.sym()) <<
"sparse problem " << it.matname()
335 <<
" (" << A.
rows() <<
"x" << A.
cols() <<
") using " <<
typeid(Solver).
name() <<
"..." << std::endl;
337 std::string
stats = solver_stats(solver);
339 std::cout <<
"INFO | " << stats << std::endl;
344 std::cout <<
"INFO | Skip sparse problem \"" << it.matname() <<
"\" (too large)" << std::endl;
371 template<
typename Solver,
typename DenseMat>
377 Index size = internal::random<int>(1,maxSize);
378 double density = (
std::max)(8./(size*size), 0.01);
381 dA.resize(size,size);
383 initSparse<Scalar>(density, dA,
A,
options);
392 template<
class Scalar>
399 template<
typename Solver>
void check_sparse_square_solving(Solver& solver,
int maxSize = 300,
int maxRealWorldSize = 100000,
bool checkDeficient =
false)
408 int rhsCols = internal::random<int>(1,16);
416 DenseVector
b = DenseVector::Random(size);
417 DenseMatrix dB(size,rhsCols);
418 SpMat
B(size,rhsCols);
419 double density = (
std::max)(8./(size*rhsCols), 0.1);
423 DenseVector dc = dB.col(0);
432 b = DenseVector::Zero(size);
436 if(checkDeficient && size>1) {
445 #ifdef TEST_REAL_CASES 449 std::string mat_folder = get_matrixfolder<Scalar>();
454 if(A.diagonal().size() <= maxRealWorldSize)
456 DenseVector
b = it.rhs();
457 DenseVector refX = it.refX();
458 std::cout <<
"INFO | Testing " << sym_to_string(it.sym()) <<
"sparse problem " << it.matname()
459 <<
" (" << A.rows() <<
"x" << A.cols() <<
") using " <<
typeid(Solver).
name() <<
"..." << std::endl;
461 std::string
stats = solver_stats(solver);
463 std::cout <<
"INFO | " << stats << std::endl;
467 std::cout <<
"INFO | SKIP sparse problem \"" << it.matname() <<
"\" (too large)" << std::endl;
488 int size = internal::random<int>(1,30);
489 dA.setRandom(size,size);
491 dA = (dA.array().abs()<0.3).select(0,dA);
492 dA.diagonal() = (dA.diagonal().array()==0).select(1,dA.diagonal());
516 template<
typename Solver,
typename DenseMat>
522 int rows = internal::random<int>(1,maxSize);
523 int cols = internal::random<int>(1,
rows);
524 double density = (
std::max)(8./(rows*cols), 0.01);
527 dA.resize(rows,cols);
529 initSparse<Scalar>(density, dA,
A,
options);
540 int rhsCols = internal::random<int>(1,16);
548 DenseVector
b = DenseVector::Random(A.rows());
549 DenseMatrix dB(A.rows(),rhsCols);
550 SpMat
B(A.rows(),rhsCols);
551 double density = (
std::max)(8./(A.rows()*rhsCols), 0.1);
561 b = DenseVector::Zero(A.rows());
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index cols() const
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
internal::traits< Derived >::Scalar Scalar
Matrix< RealScalar, Dynamic, Dynamic > M
void check_sparse_solving(Solver &solver, const typename Solver::MatrixType &A, const Rhs &b, const DenseMat &dA, const DenseRhs &db)
A matrix or vector expression mapping an existing array of data.
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Index rows() const
A base class for sparse solvers.
void check_sparse_determinant(Solver &solver, const typename Solver::MatrixType &A, const DenseMat &dA)
BiCGSTAB< SparseMatrix< double > > solver
Pose3 x2(Rot3::Ypr(0.0, 0.0, 0.0), l2)
Holds information about the various numeric (i.e. scalar) types allowed by Eigen. ...
void check_sparse_spd_determinant(Solver &solver)
bool operator()(Index, Index col, const Scalar &) const
Matrix< SCALARA, Dynamic, Dynamic > A
Matrix< SCALARB, Dynamic, Dynamic > B
void g(const string &key, int i)
Matrix< Scalar, Dynamic, 1 > DenseVector
#define VERIFY_IS_APPROX(a, b)
void check_sparse_leastsquare_solving(Solver &solver)
Base class of any sparse matrices or sparse expressions.
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
#define VERIFY_IS_EQUAL(a, b)
void check_sparse_solving_real_cases(Solver &solver, const typename Solver::MatrixType &A, const Rhs &b, const typename Solver::MatrixType &fullA, const Rhs &refX)
idx_t idx_t idx_t idx_t idx_t idx_t idx_t real_t real_t idx_t * options
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Iterator to browse matrices from a specified folder.
Matrix< Scalar, Dynamic, Dynamic > Mat
NumTraits< Scalar >::Real RealScalar
static std::stringstream ss
Eigen::SparseMatrix< double > SpMat
int generate_sparse_spd_problem(Solver &, typename Solver::MatrixType &A, typename Solver::MatrixType &halfA, DenseMat &dA, int maxSize=300)
void check_sparse_abs_determinant(Solver &solver, const typename Solver::MatrixType &A, const DenseMat &dA)
#define CALL_SUBTEST(FUNC)
void check_sparse_spd_solving(Solver &solver, int maxSize=300, int maxRealWorldSize=100000)
Index generate_sparse_square_problem(Solver &, typename Solver::MatrixType &A, DenseMat &dA, int maxSize=300, int options=ForceNonZeroDiag)
void generate_sparse_leastsquare_problem(Solver &, typename Solver::MatrixType &A, DenseMat &dA, int maxSize=300, int options=ForceNonZeroDiag)
void solve_with_guess(IterativeSolverBase< Solver > &solver, const MatrixBase< Rhs > &b, const Guess &g, Result &x)
internal::nested_eval< T, 1 >::type eval(const T &xpr)
void check_sparse_square_determinant(Solver &solver)
set noclip points set clip one set noclip two set bar set border lt lw set xdata set ydata set zdata set x2data set y2data set boxwidth set dummy x
Base class for linear iterative solvers.
NumTraits< Scalar >::Real RealScalar
Base class for all dense matrices, vectors, and expressions.
void check_sparse_square_abs_determinant(Solver &solver)
void check_sparse_square_solving(Solver &solver, int maxSize=300, int maxRealWorldSize=100000, bool checkDeficient=false)
#define EIGEN_UNUSED_VARIABLE(var)