Skip to content

improve Simplicial Cholesky analyzePattern

Reference issue

Fixes #2661 (closed)

What does this implement/fix?

Improves the runtime of the Simplicial Cholesky analyze pattern routines. Implements the bleeding edge algorithms from:

Joseph W. Liu. (1986). A compact row storage scheme for Cholesky factors using elimination trees. ACM Trans. Math. Softw. 12, 2 (June 1986), 127–148. https://doi.org/10.1145/6497.6499

and

Gilbert, J. R., Ng, E., & Peyton, B. W. (1994). An efficient algorithm to compute row and column counts for sparse Cholesky factorization. SIAM Journal on Matrix Analysis and Applications, 15(4), 1075–1091.

The benchmark problems are from https://sparse.tamu.edu/. ubuntu clang++-14 -O3 -DNDEBUG

Before:

------------------------------------------------------------------------------
Benchmark                                    Time             CPU   Iterations
------------------------------------------------------------------------------
test_sparse/crankseg_1/iterations:5 1561498721 ns   1561325246 ns            5
test_sparse/inline_1/iterations:1   7.1502e+10 ns   7.1488e+10 ns            1
test_sparse/thermal2/iterations:1   5.1971e+10 ns   5.1963e+10 ns            1
test_sparse/audikw_1/iterations:1   4.5089e+11 ns   4.5082e+11 ns            1

After:

------------------------------------------------------------------------------
Benchmark                                    Time             CPU   Iterations
------------------------------------------------------------------------------
test_sparse/crankseg_1/iterations:5   51039424 ns     51037794 ns            5
test_sparse/inline_1/iterations:1    209916927 ns    209920071 ns            1
test_sparse/thermal2/iterations:1    150026624 ns    150024657 ns            1
test_sparse/audikw_1/iterations:1    426583838 ns    426578780 ns            1

In the largest benchmark problem audikw_1, the runtime was reduced from 7.5 minutes to less than 0.5 seconds. Not too shabby!

Additional information

Google bm code:

#include <benchmark/benchmark.h>
#include <Eigen/SparseCholesky>
#include <unsupported/Eigen/SparseExtra>

using namespace Eigen;
using Scalar = double;
using StorageIndex = int;
using Mat = SparseMatrix<Scalar, ColMajor, StorageIndex>;
using Ordering = NaturalOrdering<StorageIndex>;
using Solver = SimplicialLLT<Mat, Lower, Ordering>;

static void test_sparse(benchmark::State& state, std::string name) {
  Mat A;
  loadMarket(A, name);
  for (auto s : state) {
    Solver solver;
    solver.analyzePattern(A);
    benchmark::DoNotOptimize(solver);
  }
}

BENCHMARK_CAPTURE(test_sparse, crankseg_1, std::string("crankseg_1.mtx"))->Iterations(5);
BENCHMARK_CAPTURE(test_sparse, inline_1, std::string("inline_1.mtx"))->Iterations(1);
BENCHMARK_CAPTURE(test_sparse, thermal2, std::string("thermal2.mtx"))->Iterations(1);
BENCHMARK_CAPTURE(test_sparse, audikw_1, std::string("audikw_1.mtx"))->Iterations(1);
BENCHMARK_MAIN();
Edited by Charles Schlosser

Merge request reports

Loading