Skip to content

Commit 2f4b259

Browse files
Add package version tracking and correctness checks to LinearSolveAutotune (#702)
* Add package version tracking and correctness checks to LinearSolveAutotune - Add get_package_versions() function to collect versions of all LinearSolve-related dependencies - Include versions of LinearSolve, RecursiveFactorization, CUDA, Metal, MKL_jll, blis_jll, etc. - Display package versions in telemetry output for better reproducibility - Add correctness check against LUFactorization baseline - Algorithms that fail correctness check (with lenient tolerance) are warned and excluded - Default tolerance is 1e-2 (very lenient) to catch major issues - Helps ensure benchmark results are meaningful and catch any fishy business * Add Pkg as dependency to LinearSolveAutotune Required for get_package_versions() function to use Pkg.dependencies() * Apply JuliaFormatter to modified files * Delete format_autotune.jl * Delete test_autotune_improvements.jl * Revert "Apply JuliaFormatter to modified files" This reverts commit 875dd22. --------- Co-authored-by: ChrisRackauckas <[email protected]>
1 parent aaa4260 commit 2f4b259

File tree

4 files changed

+163
-14
lines changed

4 files changed

+163
-14
lines changed

lib/LinearSolveAutotune/Project.toml

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -16,6 +16,7 @@ LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
1616
LinearSolve = "7ed4a6bd-45f5-4d41-b270-4a48e9bafcae"
1717
MKL_jll = "856f044c-d86e-5d09-b602-aeab76dc8ba7"
1818
Metal = "dde4c033-4e86-420c-a63e-0dd931031962"
19+
Pkg = "44cfe95a-1eb2-52ea-b672-e2afdf69b78f"
1920
Plots = "91a5bcdd-55d7-5caf-9e0b-520d859cae80"
2021
Preferences = "21216c6a-2e73-6563-6e65-726566657250"
2122
PrettyTables = "08abe8d2-0d0c-5749-adfa-8a2ac140af0d"
@@ -41,6 +42,7 @@ LinearAlgebra = "1"
4142
LinearSolve = "3"
4243
MKL_jll = "2025.2.0"
4344
Metal = "1"
45+
Pkg = "1"
4446
Plots = "1"
4547
Preferences = "1"
4648
PrettyTables = "2"

lib/LinearSolveAutotune/src/benchmarking.jl

Lines changed: 47 additions & 14 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,7 @@
11
# Core benchmarking functionality
22

33
using ProgressMeter
4+
using LinearAlgebra
45

56
"""
67
test_algorithm_compatibility(alg, eltype::Type, test_size::Int=4)
@@ -78,7 +79,8 @@ Benchmark the given algorithms across different matrix sizes and element types.
7879
Returns a DataFrame with results including element type information.
7980
"""
8081
function benchmark_algorithms(matrix_sizes, algorithms, alg_names, eltypes;
81-
samples = 5, seconds = 0.5, sizes = [:tiny, :small, :medium, :large])
82+
samples = 5, seconds = 0.5, sizes = [:tiny, :small, :medium, :large],
83+
check_correctness = true, correctness_tol = 1e-2)
8284

8385
# Set benchmark parameters
8486
old_params = BenchmarkTools.DEFAULT_PARAMETERS
@@ -116,6 +118,18 @@ function benchmark_algorithms(matrix_sizes, algorithms, alg_names, eltypes;
116118
A = rand(rng, eltype, n, n)
117119
b = rand(rng, eltype, n)
118120
u0 = rand(rng, eltype, n)
121+
122+
# Compute reference solution with LUFactorization if correctness check is enabled
123+
reference_solution = nothing
124+
if check_correctness
125+
try
126+
ref_prob = LinearProblem(copy(A), copy(b); u0 = copy(u0))
127+
reference_solution = solve(ref_prob, LinearSolve.LUFactorization())
128+
catch e
129+
@warn "Failed to compute reference solution with LUFactorization for size $n, eltype $eltype: $e"
130+
check_correctness = false # Disable for this size/type combination
131+
end
132+
end
119133

120134
for (alg, name) in zip(compatible_algs, compatible_names)
121135
# Update progress description
@@ -125,26 +139,45 @@ function benchmark_algorithms(matrix_sizes, algorithms, alg_names, eltypes;
125139
gflops = 0.0
126140
success = true
127141
error_msg = ""
142+
passed_correctness = true
128143

129144
try
130145
# Create the linear problem for this test
131146
prob = LinearProblem(copy(A), copy(b);
132147
u0 = copy(u0),
133148
alias = LinearAliasSpecifier(alias_A = true, alias_b = true))
134149

135-
# Warmup run
136-
solve(prob, alg)
137-
138-
# Actual benchmark
139-
bench = @benchmark solve($prob, $alg) setup=(prob = LinearProblem(
140-
copy($A), copy($b);
141-
u0 = copy($u0),
142-
alias = LinearAliasSpecifier(alias_A = true, alias_b = true)))
143-
144-
# Calculate GFLOPs
145-
min_time_sec = minimum(bench.times) / 1e9
146-
flops = luflop(n, n)
147-
gflops = flops / min_time_sec / 1e9
150+
# Warmup run and correctness check
151+
warmup_sol = solve(prob, alg)
152+
153+
# Check correctness if reference solution is available
154+
if check_correctness && reference_solution !== nothing
155+
# Compute relative error
156+
rel_error = norm(warmup_sol.u - reference_solution.u) / norm(reference_solution.u)
157+
158+
if rel_error > correctness_tol
159+
passed_correctness = false
160+
@warn "Algorithm $name failed correctness check for size $n, eltype $eltype. " *
161+
"Relative error: $(round(rel_error, sigdigits=3)) > tolerance: $correctness_tol. " *
162+
"Algorithm will be excluded from results."
163+
success = false
164+
error_msg = "Failed correctness check (rel_error = $(round(rel_error, sigdigits=3)))"
165+
end
166+
end
167+
168+
# Only benchmark if correctness check passed
169+
if passed_correctness
170+
# Actual benchmark
171+
bench = @benchmark solve($prob, $alg) setup=(prob = LinearProblem(
172+
copy($A), copy($b);
173+
u0 = copy($u0),
174+
alias = LinearAliasSpecifier(alias_A = true, alias_b = true)))
175+
176+
# Calculate GFLOPs
177+
min_time_sec = minimum(bench.times) / 1e9
178+
flops = luflop(n, n)
179+
gflops = flops / min_time_sec / 1e9
180+
end
148181

149182
catch e
150183
success = false

lib/LinearSolveAutotune/src/gpu_detection.jl

Lines changed: 99 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,7 @@
11
# GPU hardware and package detection
22

33
using CPUSummary
4+
using Pkg
45

56
"""
67
is_cuda_available()
@@ -116,9 +117,107 @@ function get_system_info()
116117
info["apple_accelerate_available"] = false
117118
end
118119

120+
# Add package versions
121+
info["package_versions"] = get_package_versions()
122+
119123
return info
120124
end
121125

126+
"""
127+
get_package_versions()
128+
129+
Get versions of LinearSolve-related packages and their dependencies.
130+
Returns a Dict with package names and versions.
131+
"""
132+
function get_package_versions()
133+
versions = Dict{String, String}()
134+
135+
# Get the current project's dependencies
136+
deps = Pkg.dependencies()
137+
138+
# List of packages we're interested in tracking
139+
important_packages = [
140+
"LinearSolve",
141+
"LinearSolveAutotune",
142+
"RecursiveFactorization",
143+
"CUDA",
144+
"Metal",
145+
"MKL_jll",
146+
"BLISBLAS",
147+
"AppleAccelerate",
148+
"SparseArrays",
149+
"KLU",
150+
"Pardiso",
151+
"MKLPardiso",
152+
"BandedMatrices",
153+
"FastLapackInterface",
154+
"HYPRE",
155+
"IterativeSolvers",
156+
"Krylov",
157+
"KrylovKit",
158+
"LinearAlgebra"
159+
]
160+
161+
# Also track JLL packages for BLAS libraries
162+
jll_packages = [
163+
"MKL_jll",
164+
"OpenBLAS_jll",
165+
"OpenBLAS32_jll",
166+
"blis_jll",
167+
"LAPACK_jll",
168+
"CompilerSupportLibraries_jll"
169+
]
170+
171+
all_packages = union(important_packages, jll_packages)
172+
173+
# Iterate through dependencies and collect versions
174+
for (uuid, dep) in deps
175+
if dep.name in all_packages
176+
if dep.version !== nothing
177+
versions[dep.name] = string(dep.version)
178+
else
179+
# Try to get version from the package itself if loaded
180+
try
181+
pkg_module = Base.loaded_modules[Base.PkgId(uuid, dep.name)]
182+
if isdefined(pkg_module, :version)
183+
versions[dep.name] = string(pkg_module.version)
184+
else
185+
versions[dep.name] = "unknown"
186+
end
187+
catch
188+
versions[dep.name] = "unknown"
189+
end
190+
end
191+
end
192+
end
193+
194+
# Try to get Julia's LinearAlgebra stdlib version
195+
try
196+
versions["LinearAlgebra"] = string(VERSION) # Stdlib version matches Julia
197+
catch
198+
versions["LinearAlgebra"] = "stdlib"
199+
end
200+
201+
# Get BLAS configuration info
202+
try
203+
blas_config = LinearAlgebra.BLAS.get_config()
204+
if hasfield(typeof(blas_config), :loaded_libs)
205+
for lib in blas_config.loaded_libs
206+
if hasfield(typeof(lib), :libname)
207+
lib_name = basename(string(lib.libname))
208+
# Extract version info if available
209+
versions["BLAS_lib"] = lib_name
210+
end
211+
end
212+
end
213+
catch
214+
# Fallback for older Julia versions
215+
versions["BLAS_vendor"] = string(LinearAlgebra.BLAS.vendor())
216+
end
217+
218+
return versions
219+
end
220+
122221
"""
123222
get_detailed_system_info()
124223

lib/LinearSolveAutotune/src/telemetry.jl

Lines changed: 15 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -172,6 +172,21 @@ function format_system_info_markdown(system_info::Dict)
172172
push!(lines, "- **CUDA Available**: $(get(system_info, "cuda_available", get(system_info, "has_cuda", false)))")
173173
# Handle both "has_metal" and "metal_available" keys
174174
push!(lines, "- **Metal Available**: $(get(system_info, "metal_available", get(system_info, "has_metal", false)))")
175+
176+
# Add package versions section
177+
if haskey(system_info, "package_versions")
178+
push!(lines, "")
179+
push!(lines, "### Package Versions")
180+
pkg_versions = system_info["package_versions"]
181+
182+
# Sort packages for consistent display
183+
sorted_packages = sort(collect(keys(pkg_versions)))
184+
185+
for pkg_name in sorted_packages
186+
version = pkg_versions[pkg_name]
187+
push!(lines, "- **$pkg_name**: $version")
188+
end
189+
end
175190

176191
return join(lines, "\n")
177192
end

0 commit comments

Comments
 (0)