LinearSolve solvers

In this page, we will benchmark several solvers provided by LinearSolve.jl for solving SteadyState and spectrum in hierarchical equations of motion approach.

using LinearSolve
using BenchmarkTools
using HierarchicalEOM
HierarchicalEOM.versioninfo()

                                   __
                                  /  \
 __     __                     __ \__/ __
|  |   |  |                   /  \    /  \
|  |   |  | ______   ______   \__/_  _\__/
|  |___|  |/  __  \ /  __  \ / '   \/     \
|   ___   |  |__)  |  /  \  |    _     _   |
|  |   |  |   ____/| (    ) |   / \   / \  |
|  |   |  |  |____ |  \__/  |  |   | |   | |
|__|   |__|\______) \______/|__|   |_|   |_|

Julia framework for Hierarchical Equations of Motion
≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡≡
Copyright © NCKU-QFORT 2022 and later.
Lead  developer : Yi-Te Huang
Other developers:
    Simon Cross, Neill Lambert, Po-Chen Kuo and Shen-Liang Yang

Package information:
====================================
HierarchicalEOM Ver. 1.6.0
LinearSolve     Ver. 2.22.1
OrdinaryDiffEq  Ver. 6.80.1
FastExpm        Ver. 1.1.0
JLD2            Ver. 0.4.48

System information:
====================================
Julia Version: 1.10.4
OS       : Linux (x86_64-linux-gnu)
CPU      : 4 × AMD EPYC 7763 64-Core Processor
Memory   : 15.606 GB
WORD_SIZE: 64
LIBM     : libopenlibm
LLVM     : libLLVM-15.0.7 (ORCJIT, znver3)
BLAS     : libopenblas64_.so (ilp64)

Here, we use the example of the single-impurity Anderson model:

ϵ = -5
U = 10
Γ = 2
μ = 0
W = 10
kT = 0.5
N = 5
tier = 2
ωlist = -10:1:10

σm = [0 1; 0 0]
σz = [1 0; 0 -1]
II = [1 0; 0 1]
d_up = kron(σm, II)
d_dn = kron(-1 * σz, σm)
Hsys = ϵ * (d_up' * d_up + d_dn' * d_dn) + U * (d_up' * d_up * d_dn' * d_dn)

bath_up = Fermion_Lorentz_Pade(d_up, Γ, μ, W, kT, N)
bath_dn = Fermion_Lorentz_Pade(d_dn, Γ, μ, W, kT, N)
bath_list = [bath_up, bath_dn]
M_even = M_Fermion(Hsys, tier, bath_list)
M_odd = M_Fermion(Hsys, tier, bath_list, ODD)
ados_s = SteadyState(M_even);
Preparing block matrices for HEOM Liouvillian superoperator (using 1 threads)...
Constructing matrix...[DONE]
Preparing block matrices for HEOM Liouvillian superoperator (using 1 threads)...
Constructing matrix...[DONE]
Solving steady state for ADOs by linear-solve method...[DONE]

LinearSolve Solver List

(click here to see the full solver list provided by LinearSolve.jl)

UMFPACKFactorization (Default solver)

This solver performs better when there is more structure to the sparsity pattern (depends on the complexity of your system and baths).

UMFPACKFactorization();

KLUFactorization

This solver performs better when there is less structure to the sparsity pattern (depends on the complexity of your system and baths).

KLUFactorization();

Julia's built-in LU factorization

LUFactorization();

A generic BICGSTAB implementation from Krylov

KrylovJL_BICGSTAB();

Pardiso

This solver is based on Intel openAPI Math Kernel Library (MKL) Pardiso

Note

Using this solver requires adding the package Pardiso.jl, i.e. using Pardiso

using Pardiso
using LinearSolve
MKLPardisoFactorize()
MKLPardisoIterate();

Solving Stationary State

Since we are using BenchmarkTools (@benchmark) in the following, we set verbose = false to disable the output message.

UMFPACKFactorization (Default solver)

@benchmark SteadyState(M_even; verbose = false)
BenchmarkTools.Trial: 746 samples with 1 evaluation.
 Range (minmax):  6.408 ms 12.067 ms   GC (min … max): 0.00% … 43.48%
 Time  (median):     6.518 ms                GC (median):    0.00%
 Time  (mean ± σ):   6.694 ms ± 444.884 μs   GC (mean ± σ):  1.47% ±  3.28%

    ▃█▆                                                       
  ▃▄████▅▄▃▂▂▂▂▁▂▂▂▁▁▁▁▁▁▂▂▂▄▄▆▅▄▃▄▂▂▂▂▁▂▂▂▁▂▁▁▁▁▁▁▁▂▂▂▂▂▃▂ ▃
  6.41 ms         Histogram: frequency by time        7.67 ms <

 Memory estimate: 9.53 MiB, allocs estimate: 150.

KLUFactorization

@benchmark SteadyState(M_even; solver = KLUFactorization(), verbose = false)
BenchmarkTools.Trial: 1421 samples with 1 evaluation.
 Range (minmax):  3.350 ms  7.715 ms   GC (min … max): 0.00% … 7.84%
 Time  (median):     3.438 ms                GC (median):    0.00%
 Time  (mean ± σ):   3.516 ms ± 260.992 μs   GC (mean ± σ):  1.37% ± 3.52%

      ▄██                                                 
  ▂▂▄▇████▆▄▃▂▂▁▁▁▁▁▁▁▂▁▁▁▁▂▁▁▂▁▂▁▂▁▁▂▁▁▁▁▁▂▁▂▂▂▂▃▃▃▃▃▃▃▃▂▂ ▃
  3.35 ms         Histogram: frequency by time        4.15 ms <

 Memory estimate: 5.85 MiB, allocs estimate: 205.

Julia's built-in generic LU factorization

@benchmark SteadyState(M_even; solver = LUFactorization(), verbose = false)
BenchmarkTools.Trial: 759 samples with 1 evaluation.
 Range (minmax):  6.354 ms  9.389 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     6.429 ms                GC (median):    0.00%
 Time  (mean ± σ):   6.584 ms ± 312.391 μs   GC (mean ± σ):  1.25% ± 2.52%

    ▇█▇                                                       
  ▅█████▇▄▄▃▃▂▂▂▂▂▁▂▁▁▁▁▁▁▁▂▁▁▁▁▁▂▂▁▁▁▁▃▃▅▆▆▅▄▄▄▃▃▂▂▁▂▁▁▁▂▂ ▃
  6.35 ms         Histogram: frequency by time         7.2 ms <

 Memory estimate: 9.40 MiB, allocs estimate: 167.

KrylovJL_BICGSTAB

@benchmark SteadyState(M_even; solver = KrylovJL_BICGSTAB(rtol = 1e-10, atol = 1e-12), verbose = false)
BenchmarkTools.Trial: 6601 samples with 1 evaluation.
 Range (minmax):  705.053 μs  3.715 ms   GC (min … max): 0.00% … 16.87%
 Time  (median):     727.105 μs                GC (median):    0.00%
 Time  (mean ± σ):   755.224 μs ± 125.125 μs   GC (mean ± σ):  2.60% ±  7.47%

  ▄█▆▅▂▁                                                 ▁ ▁▁ ▂
  ██████▇▆▄▅▃▄▁▄▄▄▃▁▁▃▆▇▇▄▅▅▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▅▇████ █
  705 μs        Histogram: log(frequency) by time       1.32 ms <

 Memory estimate: 1.85 MiB, allocs estimate: 87.

MKLPardisoFactorize

@benchmark SteadyState(M_even; solver = MKLPardisoFactorize(), verbose = false)
BenchmarkTools.Trial: 376 samples with 1 evaluation.
 Range (minmax):  11.874 ms35.927 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     13.774 ms               GC (median):    0.00%
 Time  (mean ± σ):   13.312 ms ±  1.799 ms   GC (mean ± σ):  0.20% ± 0.92%

  ▅▇▅▃▁▁                                  ▁▇▇▅▄               
  ██████▇▇▄▁▄▁▁▁▁▁▄▁▁▄▁▄▆▁▇▆▄▁▁▆█▇▄▇▄▄▆▄▄█████▇▄▁▄▁▁▁▄▁▆▁▁▄ ▇
  11.9 ms      Histogram: log(frequency) by time      14.6 ms <

 Memory estimate: 2.39 MiB, allocs estimate: 92.

MKLPardisoIterate

@benchmark SteadyState(M_even; solver = MKLPardisoIterate(), verbose = false)
BenchmarkTools.Trial: 369 samples with 1 evaluation.
 Range (minmax):  11.885 ms60.422 ms   GC (min … max): 0.00% … 0.00%
 Time  (median):     13.759 ms               GC (median):    0.00%
 Time  (mean ± σ):   13.573 ms ±  3.645 ms   GC (mean ± σ):  0.19% ± 0.88%

  ▇▄               ▅                                         
  ██▇▁▄▄▁▁▆▇▅▆▄██▄██▄▅▅▆▁▁▁▁▁▁▁▄▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▄▁▁▇ ▆
  11.9 ms      Histogram: log(frequency) by time        18 ms <

 Memory estimate: 2.39 MiB, allocs estimate: 92.

Calculate Spectrum

UMFPACKFactorization (Default solver)

@benchmark DensityOfStates(M_odd, ados_s, d_up, ωlist; verbose = false)
BenchmarkTools.Trial: 20 samples with 1 evaluation.
 Range (minmax):  214.618 ms695.153 ms   GC (min … max):  0.73% … 68.93%
 Time  (median):     219.136 ms                GC (median):     0.78%
 Time  (mean ± σ):   251.868 ms ± 106.171 ms   GC (mean ± σ):  10.62% ± 15.27%

  █     ▁▄▁▁▄▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▃ ▁
  215 ms           Histogram: frequency by time          695 ms <

 Memory estimate: 267.37 MiB, allocs estimate: 27816.

KLUFactorization

@benchmark DensityOfStates(M_odd, ados_s, d_up, ωlist; solver = KLUFactorization(), verbose = false)
BenchmarkTools.Trial: 71 samples with 1 evaluation.
 Range (minmax):  67.853 ms122.602 ms   GC (min … max): 0.77% … 0.96%
 Time  (median):     68.220 ms                GC (median):    0.87%
 Time  (mean ± σ):   70.428 ms ±   8.765 ms   GC (mean ± σ):  0.92% ± 0.18%

  █   ▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▄▁▁▁▁▁▁▁▁▁▁▁▁▁▄▁▁▁▁▁▁▁▄▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▄ ▁
  67.9 ms       Histogram: log(frequency) by time       106 ms <

 Memory estimate: 86.96 MiB, allocs estimate: 27505.

Julia's built-in LU factorization

@benchmark DensityOfStates(M_odd, ados_s, d_up, ωlist; solver = LUFactorization(), verbose = false)
BenchmarkTools.Trial: 17 samples with 1 evaluation.
 Range (minmax):  265.954 ms771.501 ms   GC (min … max):  0.93% … 65.35%
 Time  (median):     267.479 ms                GC (median):     0.84%
 Time  (mean ± σ):   307.134 ms ± 121.003 ms   GC (mean ± σ):  10.35% ± 15.65%

                                                                
  ▁▃▅▁▁▃▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▃ ▁
  266 ms           Histogram: frequency by time          772 ms <

 Memory estimate: 344.92 MiB, allocs estimate: 29300.

KrylovJL_BICGSTAB

@benchmark DensityOfStates(
    M_odd,
    ados_s,
    d_up,
    ωlist;
    solver = KrylovJL_BICGSTAB(rtol = 1e-10, atol = 1e-12),
    verbose = false,
)
BenchmarkTools.Trial: 15 samples with 1 evaluation.
 Range (minmax):  345.833 ms383.277 ms   GC (min … max): 0.17% … 0.16%
 Time  (median):     347.744 ms                GC (median):    0.17%
 Time  (mean ± σ):   350.842 ms ±   9.392 ms   GC (mean ± σ):  0.17% ± 0.07%

  ▃▃█                                                       
  ████▁▁▁▁▇▁▁▁▁▁▁▁▇▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▇ ▁
  346 ms           Histogram: frequency by time          383 ms <

 Memory estimate: 51.67 MiB, allocs estimate: 26276.

MKLPardisoFactorize

@benchmark DensityOfStates(M_odd, ados_s, d_up, ωlist; solver = MKLPardisoFactorize(), verbose = false)
BenchmarkTools.Trial: 51 samples with 1 evaluation.
 Range (minmax):  91.628 ms131.713 ms   GC (min … max): 0.00% … 0.47%
 Time  (median):     96.695 ms                GC (median):    0.60%
 Time  (mean ± σ):   99.059 ms ±   7.677 ms   GC (mean ± σ):  0.49% ± 0.23%

         ▂  █ ▄▂                                                
  ▄▁▁▄█▄██▆█████▆▄▁▁▁▁▁▄▁▁▁▁▁▁▄▁▁▁▁▁▁▁▁▁▁▄▁▁▁▁▁▁▁▁▁▄▁▁▄▁▁▁▄▄ ▁
  91.6 ms         Histogram: frequency by time          118 ms <

 Memory estimate: 43.61 MiB, allocs estimate: 25843.

MKLPardisoIterate

@benchmark DensityOfStates(M_odd, ados_s, d_up, ωlist; solver = MKLPardisoIterate(), verbose = false)
BenchmarkTools.Trial: 50 samples with 1 evaluation.
 Range (minmax):   91.291 ms147.021 ms   GC (min … max): 0.63% … 0.38%
 Time  (median):      96.391 ms                GC (median):    0.62%
 Time  (mean ± σ):   100.053 ms ±  13.049 ms   GC (mean ± σ):  0.49% ± 0.24%

     ▃▂▆                                                     
  ▄▅█████▅▅▁▄▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▄▁▄▁▄▁▁▄ ▁
  91.3 ms          Histogram: frequency by time          147 ms <

 Memory estimate: 43.61 MiB, allocs estimate: 25843.

This page was generated using Literate.jl.