11 #include <Eigen/SparseCore> 12 #include <Eigen/SparseLU> 15 template<
typename Solver,
typename Rhs,
typename Guess,
typename Result>
17 if(internal::random<bool>())
20 x = solver.
derived().solveWithGuess(b,g) + Result::Zero(x.rows(), x.cols());
25 x = solver.
derived().solveWithGuess(b.derived(),
g);
29 template<
typename Solver,
typename Rhs,
typename Guess,
typename Result>
31 if(internal::random<bool>())
32 x = solver.
derived().solve(b) + Result::Zero(x.rows(), x.cols());
37 template<
typename Solver,
typename Rhs,
typename Guess,
typename Result>
42 template<
typename Solver,
typename Rhs,
typename DenseMat,
typename DenseRhs>
47 typedef typename Mat::StorageIndex StorageIndex;
49 DenseRhs refX = dA.householderQr().solve(db);
51 Rhs x(A.cols(), b.cols());
57 std::cerr <<
"ERROR | sparse solver testing, factorization failed (" <<
typeid(Solver).
name() <<
")\n";
63 std::cerr <<
"WARNING: sparse solver testing: solving failed (" <<
typeid(Solver).
name() <<
")\n";
70 VERIFY(oldb.isApprox(b) &&
"sparse solver testing: the rhs should not be modified!");
71 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
75 VERIFY(solver.info() ==
Success &&
"solving failed when using solve_with_guess API");
76 VERIFY(oldb.isApprox(b) &&
"sparse solver testing: the rhs should not be modified!");
77 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
81 solver.analyzePattern(A);
83 VERIFY(solver.info() ==
Success &&
"factorization failed when using analyzePattern/factorize API");
85 VERIFY(solver.info() ==
Success &&
"solving failed when using analyzePattern/factorize API");
86 VERIFY(oldb.isApprox(b) &&
"sparse solver testing: the rhs should not be modified!");
87 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
93 VERIFY(solver.info() ==
Success &&
"factorization failed when using Map");
98 xm = solver.solve(bm);
99 VERIFY(solver.info() ==
Success &&
"solving failed when using Map");
100 VERIFY(oldb.isApprox(bm) &&
"sparse solver testing: the rhs should not be modified!");
101 VERIFY(xm.isApprox(refX,test_precision<Scalar>()));
109 Rhs x(b.rows(), b.cols());
112 x = solver2.solve(b);
113 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
118 DenseRhs
x(refX.rows(), refX.cols());
121 x.block(0,0,
x.rows(),
x.cols()) = solver.solve(db.block(0,0,db.rows(),db.cols()));
122 VERIFY(oldb.isApprox(db) &&
"sparse solver testing: the rhs should not be modified!");
123 VERIFY(
x.isApprox(refX,test_precision<Scalar>()));
129 A2.reserve((ArrayXf::Random(A.outerSize())+2).template cast<typename Mat::StorageIndex>().
eval());
131 Rhs x = solver.solve(b);
132 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
137 solver.compute(0.5*(A+A));
138 Rhs x = solver.solve(b);
139 VERIFY(x.isApprox(refX,test_precision<Scalar>()));
141 Solver solver2(0.5*(A+A));
142 Rhs x2 = solver2.solve(b);
143 VERIFY(x2.isApprox(refX,test_precision<Scalar>()));
149 template<
typename Scalar,
typename Rhs,
typename DenseMat,
typename DenseRhs>
153 typedef typename Mat::StorageIndex StorageIndex;
157 DenseRhs refX1 = dA.householderQr().
solve(db);
158 DenseRhs refX2 = dA.transpose().householderQr().solve(db);
159 DenseRhs refX3 = dA.adjoint().householderQr().solve(db);
171 std::cerr <<
"ERROR | sparse solver testing, factorization failed (" <<
typeid(Solver).
name() <<
")\n";
177 std::cerr <<
"WARNING | sparse solver testing: solving failed (" <<
typeid(Solver).
name() <<
")\n";
180 VERIFY(oldb.isApprox(b,0.0) &&
"sparse solver testing: the rhs should not be modified!");
181 VERIFY(x1.isApprox(refX1,test_precision<Scalar>()));
185 VERIFY(oldb.isApprox(b) &&
"sparse solver testing: the rhs should not be modified!");
186 VERIFY(x2.isApprox(refX2,test_precision<Scalar>()));
192 VERIFY(oldb.isApprox(b,0.0) &&
"sparse solver testing: the rhs should not be modified!");
193 VERIFY(x3.isApprox(refX3,test_precision<Scalar>()));
197 VERIFY(
solver.info() ==
Success &&
"solving failed when using analyzePattern/factorize API");
198 VERIFY(oldb.isApprox(b,0.0) &&
"sparse solver testing: the rhs should not be modified!");
199 VERIFY(x1.isApprox(refX1,test_precision<Scalar>()));
207 VERIFY(
solver.info() ==
Success &&
"factorization failed when using analyzePattern/factorize API");
209 x2 =
solver.transpose().solve(b);
210 x3 =
solver.adjoint().solve(b);
212 VERIFY(
solver.info() ==
Success &&
"solving failed when using analyzePattern/factorize API");
213 VERIFY(oldb.isApprox(b,0.0) &&
"sparse solver testing: the rhs should not be modified!");
214 VERIFY(x1.isApprox(refX1,test_precision<Scalar>()));
215 VERIFY(x2.isApprox(refX2,test_precision<Scalar>()));
216 VERIFY(x3.isApprox(refX3,test_precision<Scalar>()));
229 VERIFY(oldb.isApprox(bm,0.0) &&
"sparse solver testing: the rhs should not be modified!");
230 VERIFY(xm.isApprox(refX1,test_precision<Scalar>()));
238 Rhs x(b.rows(), b.cols());
241 x = solver2.solve(b);
242 VERIFY(x.isApprox(refX1,test_precision<Scalar>()));
247 DenseRhs
x(refX1.rows(), refX1.cols());
250 x.block(0,0,
x.rows(),
x.cols()) =
solver.solve(db.block(0,0,db.rows(),db.cols()));
251 VERIFY(oldb.isApprox(db,0.0) &&
"sparse solver testing: the rhs should not be modified!");
252 VERIFY(
x.isApprox(refX1,test_precision<Scalar>()));
258 A2.reserve((ArrayXf::Random(A.
outerSize())+2).template cast<typename Mat::StorageIndex>().
eval());
261 VERIFY(x.isApprox(refX1,test_precision<Scalar>()));
266 solver.compute(0.5*(A+A));
268 VERIFY(x.isApprox(refX1,test_precision<Scalar>()));
270 Solver solver2(0.5*(A+A));
271 Rhs x2 = solver2.solve(b);
272 VERIFY(x2.isApprox(refX1,test_precision<Scalar>()));
278 template<
typename Solver,
typename Rhs>
285 Rhs x(A.cols(), b.cols());
290 std::cerr <<
"ERROR | sparse solver testing, factorization failed (" <<
typeid(Solver).
name() <<
")\n";
297 std::cerr <<
"WARNING | sparse solver testing, solving failed (" <<
typeid(Solver).
name() <<
")\n";
301 RealScalar res_error = (fullA*x-
b).norm()/b.norm();
302 VERIFY( (res_error <= test_precision<Scalar>() ) &&
"sparse solver failed without noticing it");
305 if(refX.size() != 0 && (refX -
x).norm()/refX.norm() > test_precision<Scalar>())
307 std::cerr <<
"WARNING | found solution is different from the provided reference one\n";
311 template<
typename Solver,
typename DenseMat>
320 std::cerr <<
"WARNING | sparse solver testing: factorization failed (check_sparse_determinant)\n";
324 Scalar refDet = dA.determinant();
327 template<
typename Solver,
typename DenseMat>
337 std::cerr <<
"WARNING | sparse solver testing: factorization failed (check_sparse_abs_determinant)\n";
341 Scalar refDet =
abs(dA.determinant());
345 template<
typename Solver,
typename DenseMat>
352 int size = internal::random<int>(1,maxSize);
356 DenseMatrix dM(size, size);
361 dA = dM * dM.adjoint();
363 halfA.resize(size,size);
367 halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
373 #ifdef TEST_REAL_CASES 374 template<
typename Scalar>
375 inline std::string get_matrixfolder()
377 std::string mat_folder = TEST_REAL_CASES;
378 if( internal::is_same<
Scalar, std::complex<float> >::
value || internal::is_same<
Scalar, std::complex<double> >::
value )
379 mat_folder = mat_folder +
static_cast<std::string
>(
"/complex/");
381 mat_folder = mat_folder +
static_cast<std::string
>(
"/real/");
384 std::string sym_to_string(
int sym)
387 if(sym==
SPD)
return "SPD ";
390 template<
typename Derived>
393 std::stringstream
ss;
397 template<
typename Derived>
408 typedef typename Mat::StorageIndex StorageIndex;
421 int rhsCols = internal::random<int>(1,16);
423 SpMat
B(size,rhsCols);
424 DenseVector
b = DenseVector::Random(size);
425 DenseMatrix dB(size,rhsCols);
428 DenseVector dc = dB.col(0);
442 b = DenseVector::Zero(size);
448 #ifdef TEST_REAL_CASES 452 std::string mat_folder = get_matrixfolder<Scalar>();
456 if (it.sym() ==
SPD){
458 if(A.diagonal().size() <= maxRealWorldSize)
460 DenseVector
b = it.rhs();
461 DenseVector refX = it.refX();
467 halfA.template selfadjointView<Solver::UpLo>() = A.template triangularView<Eigen::Lower>().twistedBy(pnull);
469 std::cout <<
"INFO | Testing " << sym_to_string(it.sym()) <<
"sparse problem " << it.matname()
470 <<
" (" << A.
rows() <<
"x" << A.
cols() <<
") using " <<
typeid(Solver).
name() <<
"..." << std::endl;
472 std::string
stats = solver_stats(solver);
474 std::cout <<
"INFO | " << stats << std::endl;
479 std::cout <<
"INFO | Skip sparse problem \"" << it.matname() <<
"\" (too large)" << std::endl;
506 template<
typename Solver,
typename DenseMat>
512 Index size = internal::random<int>(1,maxSize);
516 dA.resize(size,size);
518 initSparse<Scalar>(density, dA,
A,
options);
527 template<
class Scalar>
534 template<
typename Solver>
void check_sparse_square_solving(Solver& solver,
int maxSize = 300,
int maxRealWorldSize = 100000,
bool checkDeficient =
false)
543 int rhsCols = internal::random<int>(1,16);
551 DenseVector
b = DenseVector::Random(size);
552 DenseMatrix dB(size,rhsCols);
553 SpMat
B(size,rhsCols);
558 DenseVector dc = dB.col(0);
570 if(checkDeficient && size>1) {
579 #ifdef TEST_REAL_CASES 583 std::string mat_folder = get_matrixfolder<Scalar>();
588 if(A.diagonal().size() <= maxRealWorldSize)
590 DenseVector
b = it.rhs();
591 DenseVector refX = it.refX();
592 std::cout <<
"INFO | Testing " << sym_to_string(it.sym()) <<
"sparse problem " << it.matname()
593 <<
" (" << A.rows() <<
"x" << A.cols() <<
") using " <<
typeid(Solver).
name() <<
"..." << std::endl;
595 std::string
stats = solver_stats(solver);
597 std::cout <<
"INFO | " << stats << std::endl;
601 std::cout <<
"INFO | SKIP sparse problem \"" << it.matname() <<
"\" (too large)" << std::endl;
622 int size = internal::random<int>(1,30);
623 dA.setRandom(size,size);
625 dA = (dA.array().abs()<0.3).select(0,dA);
626 dA.diagonal() = (dA.diagonal().array()==0).select(1,dA.diagonal());
650 template<
typename Solver,
typename DenseMat>
656 int rows = internal::random<int>(1,maxSize);
657 int cols = internal::random<int>(1,
rows);
661 dA.resize(rows,cols);
663 initSparse<Scalar>(density, dA,
A,
options);
674 int rhsCols = internal::random<int>(1,16);
682 DenseVector
b = DenseVector::Random(A.rows());
683 DenseMatrix dB(A.rows(),rhsCols);
684 SpMat
B(A.rows(),rhsCols);
695 b = DenseVector::Zero(A.rows());
Matrix< SCALARB, Dynamic, Dynamic, opt_B > B
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)
const StorageIndex * outerIndexPtr() const
A matrix or vector expression mapping an existing array of data.
A base class for sparse solvers.
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
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)
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Sparse supernodal LU factorization for general matrices.
void g(const string &key, int i)
bool operator()(Index, Index col, const Scalar &) const
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.
NumTraits< Scalar >::Real RealScalar
Pose3 x3(Rot3::Ypr(M_PI/4.0, 0.0, 0.0), l2)
static std::stringstream ss
const StorageIndex * innerIndexPtr() const
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)
#define EIGEN_TEST_MAX_SIZE
Matrix< Scalar, Dynamic, Dynamic > Mat
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
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Base class for all dense matrices, vectors, and expressions.
const Scalar * valuePtr() const
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)
void check_sparse_spd_solving(Solver &solver, int maxSize=(std::min)(300, EIGEN_TEST_MAX_SIZE), int maxRealWorldSize=100000)