sparse_solver.h
Go to the documentation of this file.
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2011 Gael Guennebaud <g.gael@free.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 #include "sparse.h"
11 #include <Eigen/SparseCore>
12 #include <Eigen/SparseLU>
13 #include <sstream>
14 
15 template <typename Solver, typename Rhs, typename Guess, typename Result>
16 void solve_with_guess(IterativeSolverBase<Solver>& solver, const MatrixBase<Rhs>& b, const Guess& g, Result& x) {
17  if (internal::random<bool>()) {
18  // With a temporary through evaluator<SolveWithGuess>
19  x = solver.derived().solveWithGuess(b, g) + Result::Zero(x.rows(), x.cols());
20  } else {
21  // direct evaluation within x through Assignment<Result,SolveWithGuess>
22  x = solver.derived().solveWithGuess(b.derived(), g);
23  }
24 }
25 
26 template <typename Solver, typename Rhs, typename Guess, typename Result>
27 void solve_with_guess(SparseSolverBase<Solver>& solver, const MatrixBase<Rhs>& b, const Guess&, Result& x) {
28  if (internal::random<bool>())
29  x = solver.derived().solve(b) + Result::Zero(x.rows(), x.cols());
30  else
31  x = solver.derived().solve(b);
32 }
33 
34 template <typename Solver, typename Rhs, typename Guess, typename Result>
36  x = solver.derived().solve(b);
37 }
38 
39 template <typename Solver, typename Rhs, typename DenseMat, typename DenseRhs>
40 void check_sparse_solving(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b, const DenseMat& dA,
41  const DenseRhs& db) {
42  typedef typename Solver::MatrixType Mat;
43  typedef typename Mat::Scalar Scalar;
44  typedef typename Mat::StorageIndex StorageIndex;
45 
46  DenseRhs refX = dA.householderQr().solve(db);
47  {
48  Rhs x(A.cols(), b.cols());
49  Rhs oldb = b;
50 
51  solver.compute(A);
52  if (solver.info() != Success) {
53  std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
54  VERIFY(solver.info() == Success);
55  }
56  x = solver.solve(b);
57  if (solver.info() != Success) {
58  std::cerr << "WARNING: sparse solver testing: solving failed (" << typeid(Solver).name() << ")\n";
59  // dump call stack:
60  g_test_level++;
61  VERIFY(solver.info() == Success);
62  g_test_level--;
63  return;
64  }
65  VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
66  VERIFY(x.isApprox(refX, test_precision<Scalar>()));
67 
68  x.setZero();
70  VERIFY(solver.info() == Success && "solving failed when using solve_with_guess API");
71  VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
72  VERIFY(x.isApprox(refX, test_precision<Scalar>()));
73 
74  x.setZero();
75  // test the analyze/factorize API
76  solver.analyzePattern(A);
77  solver.factorize(A);
78  VERIFY(solver.info() == Success && "factorization failed when using analyzePattern/factorize API");
79  x = solver.solve(b);
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>()));
83 
84  x.setZero();
85  // test with Map
87  A.rows(), A.cols(), A.nonZeros(), const_cast<StorageIndex*>(A.outerIndexPtr()),
88  const_cast<StorageIndex*>(A.innerIndexPtr()), const_cast<Scalar*>(A.valuePtr()));
89  solver.compute(Am);
90  VERIFY(solver.info() == Success && "factorization failed when using Map");
91  DenseRhs dx(refX);
92  dx.setZero();
93  Map<DenseRhs> xm(dx.data(), dx.rows(), dx.cols());
94  Map<const DenseRhs> bm(db.data(), db.rows(), db.cols());
95  xm = solver.solve(bm);
96  VERIFY(solver.info() == Success && "solving failed when using Map");
97  VERIFY(oldb.isApprox(bm) && "sparse solver testing: the rhs should not be modified!");
98  VERIFY(xm.isApprox(refX, test_precision<Scalar>()));
99 
100  // Test with a Map and non-unit stride.
101  Eigen::Matrix<Scalar, Eigen::Dynamic, Eigen::Dynamic> out(2 * xm.rows(), 2 * xm.cols());
102  out.setZero();
103  Eigen::Map<DenseRhs, 0, Stride<Eigen::Dynamic, 2>> outm(out.data(), xm.rows(), xm.cols(),
104  Stride<Eigen::Dynamic, 2>(2 * xm.rows(), 2));
105  outm = solver.solve(bm);
106  VERIFY(outm.isApprox(refX, test_precision<Scalar>()));
107  }
108 
109  // if not too large, do some extra check:
110  if (A.rows() < 2000) {
111  // test initialization ctor
112  {
113  Rhs x(b.rows(), b.cols());
114  Solver solver2(A);
115  VERIFY(solver2.info() == Success);
116  x = solver2.solve(b);
117  VERIFY(x.isApprox(refX, test_precision<Scalar>()));
118  }
119 
120  // test dense Block as the result and rhs:
121  {
122  DenseRhs x(refX.rows(), refX.cols());
123  DenseRhs oldb(db);
124  x.setZero();
125  x.block(0, 0, x.rows(), x.cols()) = solver.solve(db.block(0, 0, db.rows(), db.cols()));
126  VERIFY(oldb.isApprox(db) && "sparse solver testing: the rhs should not be modified!");
127  VERIFY(x.isApprox(refX, test_precision<Scalar>()));
128  }
129 
130  // test uncompressed inputs
131  {
132  Mat A2 = A;
133  A2.reserve((ArrayXf::Random(A.outerSize()) + 2).template cast<typename Mat::StorageIndex>().eval());
134  solver.compute(A2);
135  Rhs x = solver.solve(b);
136  VERIFY(x.isApprox(refX, test_precision<Scalar>()));
137  }
138 
139  // test expression as input
140  {
141  solver.compute(0.5 * (A + A));
142  Rhs x = solver.solve(b);
143  VERIFY(x.isApprox(refX, test_precision<Scalar>()));
144 
145  Solver solver2(0.5 * (A + A));
146  Rhs x2 = solver2.solve(b);
147  VERIFY(x2.isApprox(refX, test_precision<Scalar>()));
148  }
149  }
150 }
151 
152 // specialization of generic check_sparse_solving for SuperLU in order to also test adjoint and transpose solves
153 template <typename Scalar, typename Rhs, typename DenseMat, typename DenseRhs>
155  const typename Eigen::SparseMatrix<Scalar>& A, const Rhs& b, const DenseMat& dA,
156  const DenseRhs& db) {
157  typedef typename Eigen::SparseMatrix<Scalar> Mat;
158  typedef typename Mat::StorageIndex StorageIndex;
159  typedef typename Eigen::SparseLU<Eigen::SparseMatrix<Scalar>> Solver;
160 
161  // reference solutions computed by dense QR solver
162  DenseRhs refX1 = dA.householderQr().solve(db); // solution of A x = db
163  DenseRhs refX2 = dA.transpose().householderQr().solve(db); // solution of A^T * x = db (use transposed matrix A^T)
164  DenseRhs refX3 = dA.adjoint().householderQr().solve(db); // solution of A^* * x = db (use adjoint matrix A^*)
165 
166  {
167  Rhs x1(A.cols(), b.cols());
168  Rhs x2(A.cols(), b.cols());
169  Rhs x3(A.cols(), b.cols());
170  Rhs oldb = b;
171 
172  solver.compute(A);
173  if (solver.info() != Success) {
174  std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
175  VERIFY(solver.info() == Success);
176  }
177  x1 = solver.solve(b);
178  if (solver.info() != Success) {
179  std::cerr << "WARNING | sparse solver testing: solving failed (" << typeid(Solver).name() << ")\n";
180  return;
181  }
182  VERIFY(oldb.isApprox(b, 0.0) && "sparse solver testing: the rhs should not be modified!");
183  VERIFY(x1.isApprox(refX1, test_precision<Scalar>()));
184 
185  // test solve with transposed
186  x2 = solver.transpose().solve(b);
187  VERIFY(oldb.isApprox(b) && "sparse solver testing: the rhs should not be modified!");
188  VERIFY(x2.isApprox(refX2, test_precision<Scalar>()));
189 
190  // test solve with adjoint
191  // solver.template _solve_impl_transposed<true>(b, x3);
192  x3 = solver.adjoint().solve(b);
193  VERIFY(oldb.isApprox(b, 0.0) && "sparse solver testing: the rhs should not be modified!");
194  VERIFY(x3.isApprox(refX3, test_precision<Scalar>()));
195 
196  x1.setZero();
198  VERIFY(solver.info() == Success && "solving failed when using analyzePattern/factorize API");
199  VERIFY(oldb.isApprox(b, 0.0) && "sparse solver testing: the rhs should not be modified!");
200  VERIFY(x1.isApprox(refX1, test_precision<Scalar>()));
201 
202  x1.setZero();
203  x2.setZero();
204  x3.setZero();
205  // test the analyze/factorize API
206  solver.analyzePattern(A);
207  solver.factorize(A);
208  VERIFY(solver.info() == Success && "factorization failed when using analyzePattern/factorize API");
209  x1 = solver.solve(b);
210  x2 = solver.transpose().solve(b);
211  x3 = solver.adjoint().solve(b);
212 
213  VERIFY(solver.info() == Success && "solving failed when using analyzePattern/factorize API");
214  VERIFY(oldb.isApprox(b, 0.0) && "sparse solver testing: the rhs should not be modified!");
215  VERIFY(x1.isApprox(refX1, test_precision<Scalar>()));
216  VERIFY(x2.isApprox(refX2, test_precision<Scalar>()));
217  VERIFY(x3.isApprox(refX3, test_precision<Scalar>()));
218 
219  x1.setZero();
220  // test with Map
222  A.rows(), A.cols(), A.nonZeros(), const_cast<StorageIndex*>(A.outerIndexPtr()),
223  const_cast<StorageIndex*>(A.innerIndexPtr()), const_cast<Scalar*>(A.valuePtr()));
224  solver.compute(Am);
225  VERIFY(solver.info() == Success && "factorization failed when using Map");
226  DenseRhs dx(refX1);
227  dx.setZero();
228  Map<DenseRhs> xm(dx.data(), dx.rows(), dx.cols());
229  Map<const DenseRhs> bm(db.data(), db.rows(), db.cols());
230  xm = solver.solve(bm);
231  VERIFY(solver.info() == Success && "solving failed when using Map");
232  VERIFY(oldb.isApprox(bm, 0.0) && "sparse solver testing: the rhs should not be modified!");
233  VERIFY(xm.isApprox(refX1, test_precision<Scalar>()));
234  }
235 
236  // if not too large, do some extra check:
237  if (A.rows() < 2000) {
238  // test initialization ctor
239  {
240  Rhs x(b.rows(), b.cols());
241  Solver solver2(A);
242  VERIFY(solver2.info() == Success);
243  x = solver2.solve(b);
244  VERIFY(x.isApprox(refX1, test_precision<Scalar>()));
245  }
246 
247  // test dense Block as the result and rhs:
248  {
249  DenseRhs x(refX1.rows(), refX1.cols());
250  DenseRhs oldb(db);
251  x.setZero();
252  x.block(0, 0, x.rows(), x.cols()) = solver.solve(db.block(0, 0, db.rows(), db.cols()));
253  VERIFY(oldb.isApprox(db, 0.0) && "sparse solver testing: the rhs should not be modified!");
254  VERIFY(x.isApprox(refX1, test_precision<Scalar>()));
255  }
256 
257  // test uncompressed inputs
258  {
259  Mat A2 = A;
260  A2.reserve((ArrayXf::Random(A.outerSize()) + 2).template cast<typename Mat::StorageIndex>().eval());
261  solver.compute(A2);
262  Rhs x = solver.solve(b);
263  VERIFY(x.isApprox(refX1, test_precision<Scalar>()));
264  }
265 
266  // test expression as input
267  {
268  solver.compute(0.5 * (A + A));
269  Rhs x = solver.solve(b);
270  VERIFY(x.isApprox(refX1, test_precision<Scalar>()));
271 
272  Solver solver2(0.5 * (A + A));
273  Rhs x2 = solver2.solve(b);
274  VERIFY(x2.isApprox(refX1, test_precision<Scalar>()));
275  }
276  }
277 }
278 
279 template <typename Solver, typename Rhs>
280 void check_sparse_solving_real_cases(Solver& solver, const typename Solver::MatrixType& A, const Rhs& b,
281  const typename Solver::MatrixType& fullA, const Rhs& refX) {
282  typedef typename Solver::MatrixType Mat;
283  typedef typename Mat::Scalar Scalar;
284  typedef typename Mat::RealScalar RealScalar;
285 
286  Rhs x(A.cols(), b.cols());
287 
288  solver.compute(A);
289  if (solver.info() != Success) {
290  std::cerr << "ERROR | sparse solver testing, factorization failed (" << typeid(Solver).name() << ")\n";
291  VERIFY(solver.info() == Success);
292  }
293  x = solver.solve(b);
294 
295  if (solver.info() != Success) {
296  std::cerr << "WARNING | sparse solver testing, solving failed (" << typeid(Solver).name() << ")\n";
297  return;
298  }
299 
300  RealScalar res_error = (fullA * x - b).norm() / b.norm();
301  VERIFY((res_error <= test_precision<Scalar>()) && "sparse solver failed without noticing it");
302 
303  if (refX.size() != 0 && (refX - x).norm() / refX.norm() > test_precision<Scalar>()) {
304  std::cerr << "WARNING | found solution is different from the provided reference one\n";
305  }
306 }
307 template <typename Solver, typename DenseMat>
308 void check_sparse_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA) {
309  typedef typename Solver::MatrixType Mat;
310  typedef typename Mat::Scalar Scalar;
311 
312  solver.compute(A);
313  if (solver.info() != Success) {
314  std::cerr << "WARNING | sparse solver testing: factorization failed (check_sparse_determinant)\n";
315  return;
316  }
317 
318  Scalar refDet = dA.determinant();
319  VERIFY_IS_APPROX(refDet, solver.determinant());
320 }
321 template <typename Solver, typename DenseMat>
322 void check_sparse_abs_determinant(Solver& solver, const typename Solver::MatrixType& A, const DenseMat& dA) {
323  using std::abs;
324  typedef typename Solver::MatrixType Mat;
325  typedef typename Mat::Scalar Scalar;
326 
327  solver.compute(A);
328  if (solver.info() != Success) {
329  std::cerr << "WARNING | sparse solver testing: factorization failed (check_sparse_abs_determinant)\n";
330  return;
331  }
332 
333  Scalar refDet = abs(dA.determinant());
334  VERIFY_IS_APPROX(refDet, solver.absDeterminant());
335 }
336 
337 template <typename Solver, typename DenseMat>
338 int generate_sparse_spd_problem(Solver&, typename Solver::MatrixType& A, typename Solver::MatrixType& halfA,
339  DenseMat& dA, int maxSize = 300) {
340  typedef typename Solver::MatrixType Mat;
341  typedef typename Mat::Scalar Scalar;
343 
344  int size = internal::random<int>(1, maxSize);
345  double density = (std::max)(8. / static_cast<double>(size * size), 0.01);
346 
347  Mat M(size, size);
348  DenseMatrix dM(size, size);
349 
350  initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
351 
352  A = M * M.adjoint();
353  dA = dM * dM.adjoint();
354 
355  halfA.resize(size, size);
356  if (Solver::UpLo == (Lower | Upper))
357  halfA = A;
358  else
359  halfA.template selfadjointView<Solver::UpLo>().rankUpdate(M);
360 
361  return size;
362 }
363 
364 #ifdef TEST_REAL_CASES
365 template <typename Scalar>
366 inline std::string get_matrixfolder() {
367  std::string mat_folder = TEST_REAL_CASES;
368  if (internal::is_same<Scalar, std::complex<float>>::value || internal::is_same<Scalar, std::complex<double>>::value)
369  mat_folder = mat_folder + static_cast<std::string>("/complex/");
370  else
371  mat_folder = mat_folder + static_cast<std::string>("/real/");
372  return mat_folder;
373 }
374 std::string sym_to_string(int sym) {
375  if (sym == Symmetric) return "Symmetric ";
376  if (sym == SPD) return "SPD ";
377  return "";
378 }
379 template <typename Derived>
380 std::string solver_stats(const IterativeSolverBase<Derived>& solver) {
381  std::stringstream ss;
382  ss << solver.iterations() << " iters, error: " << solver.error();
383  return ss.str();
384 }
385 template <typename Derived>
386 std::string solver_stats(const SparseSolverBase<Derived>& /*solver*/) {
387  return "";
388 }
389 #endif
390 
391 template <typename Solver>
392 void check_sparse_spd_solving(Solver& solver, int maxSize = (std::min)(300, EIGEN_TEST_MAX_SIZE),
393  int maxRealWorldSize = 100000) {
394  typedef typename Solver::MatrixType Mat;
395  typedef typename Mat::Scalar Scalar;
396  typedef typename Mat::StorageIndex StorageIndex;
401 
402  // generate the problem
403  Mat A, halfA;
404  DenseMatrix dA;
405  for (int i = 0; i < g_repeat; i++) {
406  int size = generate_sparse_spd_problem(solver, A, halfA, dA, maxSize);
407 
408  // generate the right hand sides
409  int rhsCols = internal::random<int>(1, 16);
410  double density = (std::max)(8. / static_cast<double>(size * rhsCols), 0.1);
411  SpMat B(size, rhsCols);
412  DenseVector b = DenseVector::Random(size);
413  DenseMatrix dB(size, rhsCols);
414  initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
415  SpVec c = B.col(0);
416  DenseVector dc = dB.col(0);
417 
419  CALL_SUBTEST(check_sparse_solving(solver, halfA, b, dA, b));
420  CALL_SUBTEST(check_sparse_solving(solver, A, dB, dA, dB));
421  CALL_SUBTEST(check_sparse_solving(solver, halfA, dB, dA, dB));
423  CALL_SUBTEST(check_sparse_solving(solver, halfA, B, dA, dB));
425  CALL_SUBTEST(check_sparse_solving(solver, halfA, c, dA, dc));
426 
427  // check only once
428  if (i == 0) {
430  check_sparse_solving(solver, A, b, dA, b);
431  }
432  }
433 
434  // First, get the folder
435 #ifdef TEST_REAL_CASES
436  // Test real problems with double precision only
437  if (internal::is_same<typename NumTraits<Scalar>::Real, double>::value) {
438  std::string mat_folder = get_matrixfolder<Scalar>();
439  MatrixMarketIterator<Scalar> it(mat_folder);
440  for (; it; ++it) {
441  if (it.sym() == SPD) {
442  A = it.matrix();
443  if (A.diagonal().size() <= maxRealWorldSize) {
444  DenseVector b = it.rhs();
445  DenseVector refX = it.refX();
447  halfA.resize(A.rows(), A.cols());
448  if (Solver::UpLo == (Lower | Upper))
449  halfA = A;
450  else
451  halfA.template selfadjointView<Solver::UpLo>() = A.template triangularView<Eigen::Lower>().twistedBy(pnull);
452 
453  std::cout << "INFO | Testing " << sym_to_string(it.sym()) << "sparse problem " << it.matname() << " ("
454  << A.rows() << "x" << A.cols() << ") using " << typeid(Solver).name() << "..." << std::endl;
456  std::string stats = solver_stats(solver);
457  if (stats.size() > 0) std::cout << "INFO | " << stats << std::endl;
459  } else {
460  std::cout << "INFO | Skip sparse problem \"" << it.matname() << "\" (too large)" << std::endl;
461  }
462  }
463  }
464  }
465 #else
466  EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
467 #endif
468 }
469 
470 template <typename Solver>
472  typedef typename Solver::MatrixType Mat;
473  typedef typename Mat::Scalar Scalar;
475 
476  // generate the problem
477  Mat A, halfA;
478  DenseMatrix dA;
479  generate_sparse_spd_problem(solver, A, halfA, dA, 30);
480 
481  for (int i = 0; i < g_repeat; i++) {
483  check_sparse_determinant(solver, halfA, dA);
484  }
485 }
486 
487 template <typename Solver, typename DenseMat>
489  DenseMat& dA, int maxSize = 300) {
490  typedef typename Solver::MatrixType Mat;
491  typedef typename Mat::Scalar Scalar;
493 
494  int size = internal::random<int>(1, maxSize);
495  double density = (std::max)(8. / static_cast<double>(size * size), 0.01);
496 
497  Mat M(size, size);
498  DenseMatrix dM(size, size);
499 
500  initSparse<Scalar>(density, dM, M, ForceNonZeroDiag);
501 
502  A = M * M.transpose();
503  dA = dM * dM.transpose();
504 
505  halfA.resize(size, size);
506  if (Solver::UpLo == (Lower | Upper))
507  halfA = A;
508  else
509  halfA = A.template triangularView<Solver::UpLo>();
510 
511  return size;
512 }
513 
514 template <typename Solver>
516  int maxRealWorldSize = 100000) {
517  typedef typename Solver::MatrixType Mat;
518  typedef typename Mat::Scalar Scalar;
519  typedef typename Mat::StorageIndex StorageIndex;
524 
525  // generate the problem
526  Mat A, halfA;
527  DenseMatrix dA;
528  for (int i = 0; i < g_repeat; i++) {
529  int size = generate_sparse_nonhermitian_problem(solver, A, halfA, dA, maxSize);
530 
531  // generate the right hand sides
532  int rhsCols = internal::random<int>(1, 16);
533  double density = (std::max)(8. / static_cast<double>(size * rhsCols), 0.1);
534  SpMat B(size, rhsCols);
535  DenseVector b = DenseVector::Random(size);
536  DenseMatrix dB(size, rhsCols);
537  initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
538  SpVec c = B.col(0);
539  DenseVector dc = dB.col(0);
540 
542  CALL_SUBTEST(check_sparse_solving(solver, halfA, b, dA, b));
543  CALL_SUBTEST(check_sparse_solving(solver, A, dB, dA, dB));
544  CALL_SUBTEST(check_sparse_solving(solver, halfA, dB, dA, dB));
546  CALL_SUBTEST(check_sparse_solving(solver, halfA, B, dA, dB));
548  CALL_SUBTEST(check_sparse_solving(solver, halfA, c, dA, dc));
549 
550  // check only once
551  if (i == 0) {
553  check_sparse_solving(solver, A, b, dA, b);
554  }
555  }
556 
557  EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
558 }
559 
560 template <typename Solver>
562  typedef typename Solver::MatrixType Mat;
563  typedef typename Mat::Scalar Scalar;
565 
566  // generate the problem
567  Mat A, halfA;
568  DenseMatrix dA;
570 
571  for (int i = 0; i < g_repeat; i++) {
573  check_sparse_determinant(solver, halfA, dA);
574  }
575 }
576 
577 template <typename Solver>
579  typedef typename Solver::MatrixType Mat;
580 
581  Mat A(1, 1);
582  solver.compute(A);
584 }
585 
586 template <typename Solver, typename DenseMat>
587 Index generate_sparse_square_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300,
588  int options = ForceNonZeroDiag) {
589  typedef typename Solver::MatrixType Mat;
590  typedef typename Mat::Scalar Scalar;
591 
592  Index size = internal::random<int>(1, maxSize);
593  double density = (std::max)(8. / static_cast<double>(size * size), 0.01);
594 
595  A.resize(size, size);
596  dA.resize(size, size);
597 
598  initSparse<Scalar>(density, dA, A, options);
599 
600  return size;
601 }
602 
603 struct prune_column {
606  template <class Scalar>
607  bool operator()(Index, Index col, const Scalar&) const {
608  return col != m_col;
609  }
610 };
611 
612 template <typename Solver>
613 void check_sparse_square_solving(Solver& solver, int maxSize = 300, int maxRealWorldSize = 100000,
614  bool checkDeficient = false) {
615  typedef typename Solver::MatrixType Mat;
616  typedef typename Mat::Scalar Scalar;
621 
622  int rhsCols = internal::random<int>(1, 16);
623 
624  Mat A;
625  DenseMatrix dA;
626  for (int i = 0; i < g_repeat; i++) {
628 
629  A.makeCompressed();
630  DenseVector b = DenseVector::Random(size);
631  DenseMatrix dB(size, rhsCols);
632  SpMat B(size, rhsCols);
633  double density = (std::max)(8. / double(size * rhsCols), 0.1);
634  initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
635  B.makeCompressed();
636  SpVec c = B.col(0);
637  DenseVector dc = dB.col(0);
639  CALL_SUBTEST(check_sparse_solving(solver, A, dB, dA, dB));
642 
643  // check only once
644  if (i == 0) {
646  }
647  // regression test for Bug 792 (structurally rank deficient matrices):
648  if (checkDeficient && size > 1) {
649  Index col = internal::random<int>(0, int(size - 1));
650  A.prune(prune_column(col));
651  solver.compute(A);
653  }
654  }
655 
656  // First, get the folder
657 #ifdef TEST_REAL_CASES
658  // Test real problems with double precision only
659  if (internal::is_same<typename NumTraits<Scalar>::Real, double>::value) {
660  std::string mat_folder = get_matrixfolder<Scalar>();
661  MatrixMarketIterator<Scalar> it(mat_folder);
662  for (; it; ++it) {
663  A = it.matrix();
664  if (A.diagonal().size() <= maxRealWorldSize) {
665  DenseVector b = it.rhs();
666  DenseVector refX = it.refX();
667  std::cout << "INFO | Testing " << sym_to_string(it.sym()) << "sparse problem " << it.matname() << " ("
668  << A.rows() << "x" << A.cols() << ") using " << typeid(Solver).name() << "..." << std::endl;
670  std::string stats = solver_stats(solver);
671  if (stats.size() > 0) std::cout << "INFO | " << stats << std::endl;
672  } else {
673  std::cout << "INFO | SKIP sparse problem \"" << it.matname() << "\" (too large)" << std::endl;
674  }
675  }
676  }
677 #else
678  EIGEN_UNUSED_VARIABLE(maxRealWorldSize);
679 #endif
680 }
681 
682 template <typename Solver>
684  typedef typename Solver::MatrixType Mat;
685  typedef typename Mat::Scalar Scalar;
687 
688  for (int i = 0; i < g_repeat; i++) {
689  // generate the problem
690  Mat A;
691  DenseMatrix dA;
692 
693  int size = internal::random<int>(1, 30);
694  dA.setRandom(size, size);
695 
696  dA = (dA.array().abs() < 0.3).select(0, dA);
697  dA.diagonal() = (dA.diagonal().array() == 0).select(1, dA.diagonal());
698  A = dA.sparseView();
699  A.makeCompressed();
700 
702  }
703 }
704 
705 template <typename Solver>
707  typedef typename Solver::MatrixType Mat;
708  typedef typename Mat::Scalar Scalar;
710 
711  for (int i = 0; i < g_repeat; i++) {
712  // generate the problem
713  Mat A;
714  DenseMatrix dA;
716  A.makeCompressed();
718  }
719 }
720 
721 template <typename Solver, typename DenseMat>
722 void generate_sparse_leastsquare_problem(Solver&, typename Solver::MatrixType& A, DenseMat& dA, int maxSize = 300,
723  int options = ForceNonZeroDiag) {
724  typedef typename Solver::MatrixType Mat;
725  typedef typename Mat::Scalar Scalar;
726 
727  int rows = internal::random<int>(1, maxSize);
728  int cols = internal::random<int>(1, rows);
729  double density = (std::max)(8. / (rows * cols), 0.01);
730 
731  A.resize(rows, cols);
732  dA.resize(rows, cols);
733 
734  initSparse<Scalar>(density, dA, A, options);
735 }
736 
737 template <typename Solver>
739  typedef typename Solver::MatrixType Mat;
740  typedef typename Mat::Scalar Scalar;
744 
745  int rhsCols = internal::random<int>(1, 16);
746 
747  Mat A;
748  DenseMatrix dA;
749  for (int i = 0; i < g_repeat; i++) {
751 
752  A.makeCompressed();
753  DenseVector b = DenseVector::Random(A.rows());
754  DenseMatrix dB(A.rows(), rhsCols);
755  SpMat B(A.rows(), rhsCols);
756  double density = (std::max)(8. / (A.rows() * rhsCols), 0.1);
757  initSparse<Scalar>(density, dB, B, ForceNonZeroDiag);
758  B.makeCompressed();
759  check_sparse_solving(solver, A, b, dA, b);
760  check_sparse_solving(solver, A, dB, dA, dB);
761  check_sparse_solving(solver, A, B, dA, dB);
762 
763  // check only once
764  if (i == 0) {
765  b = DenseVector::Zero(A.rows());
766  check_sparse_solving(solver, A, b, dA, b);
767  }
768  }
769 }
AnnoyingScalar abs(const AnnoyingScalar &x)
Definition: AnnoyingScalar.h:135
Matrix< Scalar, Dynamic, Dynamic > DenseMatrix
Definition: BenchSparseUtil.h:23
Matrix< Scalar, Dynamic, 1 > DenseVector
Definition: BenchSparseUtil.h:24
BiCGSTAB< SparseMatrix< double > > solver
Definition: BiCGSTAB_simple.cpp:5
int i
Definition: BiCGSTAB_step_by_step.cpp:9
#define EIGEN_UNUSED_VARIABLE(var)
Definition: Macros.h:966
m col(1)
int rows
Definition: Tutorial_commainit_02.cpp:1
int cols
Definition: Tutorial_commainit_02.cpp:1
Eigen::SparseMatrix< double > SpMat
Definition: Tutorial_sparse_example.cpp:5
Scalar Scalar int size
Definition: benchVecAdd.cpp:17
Scalar * b
Definition: benchVecAdd.cpp:17
SCALAR Scalar
Definition: bench_gemm.cpp:45
Matrix< SCALARA, Dynamic, Dynamic, opt_A > A
Definition: bench_gemm.cpp:47
Matrix< RealScalar, Dynamic, Dynamic > M
Definition: bench_gemm.cpp:50
NumTraits< Scalar >::Real RealScalar
Definition: bench_gemm.cpp:46
Matrix< SCALARB, Dynamic, Dynamic, opt_B > B
Definition: bench_gemm.cpp:48
MatrixXf MatrixType
Definition: benchmark-blocking-sizes.cpp:52
#define EIGEN_TEST_MAX_SIZE
Definition: boostmultiprec.cpp:16
Base class for linear iterative solvers.
Definition: IterativeSolverBase.h:124
A matrix or vector expression mapping an existing array of data.
Definition: Map.h:96
Base class for all dense matrices, vectors, and expressions.
Definition: MatrixBase.h:52
Iterator to browse matrices from a specified folder.
Definition: MatrixMarketIterator.h:42
VectorType & refX()
Definition: MatrixMarketIterator.h:132
std::string & matname()
Definition: MatrixMarketIterator.h:147
VectorType & rhs()
Definition: MatrixMarketIterator.h:103
MatrixType & matrix()
Definition: MatrixMarketIterator.h:70
int sym()
Definition: MatrixMarketIterator.h:149
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
NumTraits< Scalar >::Real RealScalar
Definition: PlainObjectBase.h:130
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index cols() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:192
internal::traits< Derived >::Scalar Scalar
Definition: PlainObjectBase.h:127
Derived & setRandom(Index size)
Definition: Random.h:147
EIGEN_DEVICE_FUNC constexpr EIGEN_STRONG_INLINE void resize(Index rows, Index cols)
Definition: PlainObjectBase.h:294
EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE EIGEN_CONSTEXPR Index rows() const EIGEN_NOEXCEPT
Definition: PlainObjectBase.h:191
Sparse supernodal LU factorization for general matrices.
Definition: SparseLU.h:151
Base class of any sparse matrices or sparse expressions.
Definition: SparseMatrixBase.h:30
A base class for sparse solvers.
Definition: SparseSolverBase.h:67
const Solve< Derived, Rhs > solve(const MatrixBase< Rhs > &b) const
Definition: SparseSolverBase.h:84
a sparse vector class
Definition: SparseVector.h:62
Holds strides information for Map.
Definition: Stride.h:55
Definition: matrices.h:74
#define min(a, b)
Definition: datatypes.h:22
#define max(a, b)
Definition: datatypes.h:23
Matrix< Scalar, Dynamic, Dynamic > Mat
Definition: gemm_common.h:15
@ Symmetric
Definition: Constants.h:229
@ Lower
Definition: Constants.h:211
@ Upper
Definition: Constants.h:213
@ NumericalIssue
Definition: Constants.h:442
@ Success
Definition: Constants.h:440
#define VERIFY_IS_APPROX(a, b)
Definition: integer_types.cpp:13
#define VERIFY(a)
Definition: main.h:362
#define CALL_SUBTEST(FUNC)
Definition: main.h:382
#define VERIFY_IS_EQUAL(a, b)
Definition: main.h:367
@ Rhs
Definition: TensorContractionMapper.h:20
squared absolute value
Definition: GlobalFunctions.h:87
@ SPD
Definition: MatrixMarketIterator.h:19
static int g_repeat
Definition: main.h:191
EIGEN_DEFAULT_DENSE_INDEX_TYPE Index
The Index type as used for the API.
Definition: Meta.h:83
static int g_test_level
Definition: main.h:190
Vector< double > x1(const Vector< double > &coord)
Cartesian coordinates centered at the point (0.5,1)
Definition: poisson/poisson_with_singularity/two_d_poisson.cc:86
Vector< double > x2(const Vector< double > &coord)
Cartesian coordinates centered at the point (1.5,1)
Definition: poisson/poisson_with_singularity/two_d_poisson.cc:102
density
Definition: UniformPSDSelfTest.py:19
int c
Definition: calibrate.py:100
std::string string(const unsigned &i)
Definition: oomph_definitions.cc:286
double Zero
Definition: pseudosolid_node_update_elements.cc:35
list x
Definition: plotDoE.py:28
string name
Definition: plotDoE.py:33
@ ForceNonZeroDiag
Definition: sparse.h:32
int generate_sparse_spd_problem(Solver &, typename Solver::MatrixType &A, typename Solver::MatrixType &halfA, DenseMat &dA, int maxSize=300)
Definition: sparse_solver.h:338
void check_sparse_square_solving(Solver &solver, int maxSize=300, int maxRealWorldSize=100000, bool checkDeficient=false)
Definition: sparse_solver.h:613
void generate_sparse_leastsquare_problem(Solver &, typename Solver::MatrixType &A, DenseMat &dA, int maxSize=300, int options=ForceNonZeroDiag)
Definition: sparse_solver.h:722
void check_sparse_nonhermitian_solving(Solver &solver, int maxSize=(std::min)(300, EIGEN_TEST_MAX_SIZE), int maxRealWorldSize=100000)
Definition: sparse_solver.h:515
void check_sparse_zero_matrix(Solver &solver)
Definition: sparse_solver.h:578
Index generate_sparse_square_problem(Solver &, typename Solver::MatrixType &A, DenseMat &dA, int maxSize=300, int options=ForceNonZeroDiag)
Definition: sparse_solver.h:587
void check_sparse_solving(Solver &solver, const typename Solver::MatrixType &A, const Rhs &b, const DenseMat &dA, const DenseRhs &db)
Definition: sparse_solver.h:40
void check_sparse_square_abs_determinant(Solver &solver)
Definition: sparse_solver.h:706
void check_sparse_spd_determinant(Solver &solver)
Definition: sparse_solver.h:471
int generate_sparse_nonhermitian_problem(Solver &, typename Solver::MatrixType &A, typename Solver::MatrixType &halfA, DenseMat &dA, int maxSize=300)
Definition: sparse_solver.h:488
void check_sparse_nonhermitian_determinant(Solver &solver)
Definition: sparse_solver.h:561
void solve_with_guess(IterativeSolverBase< Solver > &solver, const MatrixBase< Rhs > &b, const Guess &g, Result &x)
Definition: sparse_solver.h:16
void check_sparse_spd_solving(Solver &solver, int maxSize=(std::min)(300, EIGEN_TEST_MAX_SIZE), int maxRealWorldSize=100000)
Definition: sparse_solver.h:392
void check_sparse_determinant(Solver &solver, const typename Solver::MatrixType &A, const DenseMat &dA)
Definition: sparse_solver.h:308
void check_sparse_abs_determinant(Solver &solver, const typename Solver::MatrixType &A, const DenseMat &dA)
Definition: sparse_solver.h:322
void check_sparse_leastsquare_solving(Solver &solver)
Definition: sparse_solver.h:738
void check_sparse_solving_real_cases(Solver &solver, const typename Solver::MatrixType &A, const Rhs &b, const typename Solver::MatrixType &fullA, const Rhs &refX)
Definition: sparse_solver.h:280
void check_sparse_square_determinant(Solver &solver)
Definition: sparse_solver.h:683
Holds information about the various numeric (i.e. scalar) types allowed by Eigen.
Definition: NumTraits.h:217
Definition: sparse_solver.h:603
Index m_col
Definition: sparse_solver.h:604
bool operator()(Index, Index col, const Scalar &) const
Definition: sparse_solver.h:607
prune_column(Index col)
Definition: sparse_solver.h:605
std::ofstream out("Result.txt")