Skip to content

Commit 29e13c4

Browse files
Add a standard LU performance benchmark
Same as SciML/SciMLBenchmarks.jl#623
1 parent aabe2f2 commit 29e13c4

File tree

1 file changed

+51
-0
lines changed

1 file changed

+51
-0
lines changed

benchmarks/lu.jl

Lines changed: 51 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,51 @@
1+
using BenchmarkTools, Random, VectorizationBase
2+
using LinearAlgebra, LinearSolve, MKL_jll
3+
nc = min(Int(VectorizationBase.num_cores()), Threads.nthreads())
4+
BLAS.set_num_threads(nc)
5+
BenchmarkTools.DEFAULT_PARAMETERS.seconds = 0.5
6+
7+
function luflop(m, n = m; innerflop = 2)
8+
sum(1:min(m, n)) do k
9+
invflop = 1
10+
scaleflop = isempty((k + 1):m) ? 0 : sum((k + 1):m)
11+
updateflop = isempty((k + 1):n) ? 0 :
12+
sum((k + 1):n) do j
13+
isempty((k + 1):m) ? 0 : sum((k + 1):m) do i
14+
innerflop
15+
end
16+
end
17+
invflop + scaleflop + updateflop
18+
end
19+
end
20+
21+
algs = [LUFactorization(), GenericLUFactorization(), RFLUFactorization(), MKLLUFactorization(), FastLUFactorization(), SimpleLUFactorization()]
22+
res = [Float64[] for i in 1:length(algs)]
23+
24+
ns = 4:8:500
25+
for i in 1:length(ns)
26+
n = ns[i]
27+
@info "$n × $n"
28+
rng = MersenneTwister(123)
29+
global A = rand(rng, n, n)
30+
global b = rand(rng, n)
31+
global u0= rand(rng, n)
32+
33+
for j in 1:length(algs)
34+
bt = @belapsed solve(prob, $(algs[j])).u setup=(prob = LinearProblem(copy(A), copy(b); u0 = copy(u0), alias_A=true, alias_b=true))
35+
push!(res[j], luflop(n) / bt / 1e9)
36+
end
37+
end
38+
39+
using Plots
40+
__parameterless_type(T) = Base.typename(T).wrapper
41+
parameterless_type(x) = __parameterless_type(typeof(x))
42+
parameterless_type(::Type{T}) where {T} = __parameterless_type(T)
43+
44+
p = plot(ns, res[1]; ylabel = "GFLOPs", xlabel = "N", title = "GFLOPs for NxN LU Factorization", label = string(Symbol(parameterless_type(algs[1]))), legend=:outertopright)
45+
for i in 2:length(res)
46+
plot!(p, ns, res[i]; label = string(Symbol(parameterless_type(algs[i]))))
47+
end
48+
p
49+
50+
savefig("lubench.png")
51+
savefig("lubench.pdf")

0 commit comments

Comments
 (0)