Skip to content

Commit c9549f2

Browse files
ChrisRackauckas-ClaudeChrisRackauckasclaude
authored
Add @static if Base.USE_GPL_LIBS guards for SuiteSparse usage (#797)
* Add @static if Base.USE_GPL_LIBS guards for SuiteSparse usage This commit adds conditional compilation guards to make LinearSolve.jl GPL-safe when Julia is built without GPL libraries (USE_GPL_LIBS=0). Changes: - Guard UMFPACKFactorization struct and dispatches - Guard CHOLMODFactorization struct and dispatches - Guard KLUFactorization dispatches (KLU is GPL) - Guard SPQR (sparse QR) references in _ldiv! methods - Add fallback error messages when GPL algorithms are used without GPL libs - Conditionally export GPL-dependent types - Add alternative defaultalg for sparse matrices without GPL libs When USE_GPL_LIBS=0: - Sparse LU falls back to SparspakFactorization (MIT licensed) - Symmetric sparse uses CholeskyFactorization instead of CHOLMOD - Sparse QR uses base Julia QR without SPQR 🤖 Generated with [Claude Code](https://claude.com/claude-code) Co-Authored-By: Claude <noreply@anthropic.com> * Refine GPL guards: keep init_cacheval stubs and exports Changes: - Move `nothing`-returning init_cacheval methods outside GPL guards so they exist as fallbacks even without GPL libs - Remove unnecessary QR/QRCompactWY fallbacks in _ldiv! - these are in base Julia, not GPL - Keep SPQR.QRSparse and CHOLMOD.Factor _ldiv! methods guarded - Remove conditional exports - algorithms can be exported but will error nicely when used without GPL libs - More conservative guards: only guard code that directly uses SparseArrays.UMFPACK, SparseArrays.CHOLMOD, etc. This ensures better compatibility and clearer error messages. 🤖 Generated with [Claude Code](https://claude.com/claude-code) Co-Authored-By: Claude <noreply@anthropic.com> * Fix GPL guards: revert factorization.jl, add solve! error branches Changes: - Revert all changes to src/factorization.jl - algorithms are always defined - Add `else` branches in extension with informative solve! errors for GPL algorithms - Fix `handle_sparsematrixcsc_lu` to use inline @static if instead of branches - Guard KLUFactorization solve! and KLU-specific init_cacheval methods - Keep fallback init_cacheval methods (returning `nothing`) outside guards When USE_GPL_LIBS=false: - Algorithms are still defined and exported - Using them gives clear error: "requires GPL libraries, rebuild Julia with USE_GPL_LIBS=1" - init_cacheval fallbacks exist but specific GPL type methods are unavailable 🤖 Generated with [Claude Code](https://claude.com/claude-code) Co-Authored-By: Claude <noreply@anthropic.com> --------- Co-authored-by: ChrisRackauckas <me@chrisrackauckas.com> Co-authored-by: Claude <noreply@anthropic.com>
1 parent dc4d9db commit c9549f2

File tree

2 files changed

+131
-41
lines changed

2 files changed

+131
-41
lines changed

ext/LinearSolveSparseArraysExt.jl

Lines changed: 116 additions & 38 deletions
Original file line numberDiff line numberDiff line change
@@ -39,14 +39,25 @@ function LinearSolve.init_cacheval(alg::RFLUFactorization,
3939
end
4040

4141
function LinearSolve.handle_sparsematrixcsc_lu(A::AbstractSparseMatrixCSC)
42-
lu(SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A), nonzeros(A)),
43-
check = false)
42+
@static if Base.USE_GPL_LIBS
43+
lu(SparseMatrixCSC(size(A)..., getcolptr(A), rowvals(A), nonzeros(A)),
44+
check = false)
45+
else
46+
error("Sparse LU factorization requires GPL libraries (UMFPACK). Use `using Sparspak` for a non-GPL alternative or rebuild Julia with USE_GPL_LIBS=1")
47+
end
4448
end
4549

50+
@static if Base.USE_GPL_LIBS
4651
function LinearSolve.defaultalg(
4752
A::Symmetric{<:BLASELTYPES, <:SparseMatrixCSC}, b, ::OperatorAssumptions{Bool})
4853
LinearSolve.DefaultLinearSolver(LinearSolve.DefaultAlgorithmChoice.CHOLMODFactorization)
4954
end
55+
else
56+
function LinearSolve.defaultalg(
57+
A::Symmetric{<:BLASELTYPES, <:SparseMatrixCSC}, b, ::OperatorAssumptions{Bool})
58+
LinearSolve.DefaultLinearSolver(LinearSolve.DefaultAlgorithmChoice.CholeskyFactorization)
59+
end
60+
end # @static if Base.USE_GPL_LIBS
5061

5162
function LinearSolve.defaultalg(A::AbstractSparseMatrixCSC{Tv, Ti}, b,
5263
assump::OperatorAssumptions{Bool}) where {Tv, Ti}
@@ -71,9 +82,13 @@ function LinearSolve.init_cacheval(alg::GenericFactorization,
7182
LinearSolve.do_factorization(alg, newA, b, u)
7283
end
7384

85+
@static if Base.USE_GPL_LIBS
86+
7487
const PREALLOCATED_UMFPACK = SparseArrays.UMFPACK.UmfpackLU(SparseMatrixCSC(0, 0, [1],
7588
Int[], Float64[]))
7689

90+
end # @static if Base.USE_GPL_LIBS
91+
7792
function LinearSolve.init_cacheval(
7893
alg::LUFactorization, A::AbstractSparseArray{<:Number, <:Integer}, b, u,
7994
Pl, Pr,
@@ -98,6 +113,8 @@ function LinearSolve.init_cacheval(
98113
nothing
99114
end
100115

116+
@static if Base.USE_GPL_LIBS
117+
101118
function LinearSolve.init_cacheval(
102119
alg::LUFactorization, A::AbstractSparseArray{Float64, Int64}, b, u,
103120
Pl, Pr,
@@ -132,6 +149,8 @@ function LinearSolve.init_cacheval(
132149
end
133150
end
134151

152+
end # @static if Base.USE_GPL_LIBS
153+
135154
function LinearSolve.init_cacheval(
136155
alg::LUFactorization, A::LinearSolve.GPUArraysCore.AnyGPUArray, b, u,
137156
Pl, Pr,
@@ -140,14 +159,6 @@ function LinearSolve.init_cacheval(
140159
ArrayInterface.lu_instance(A)
141160
end
142161

143-
function LinearSolve.init_cacheval(
144-
alg::UMFPACKFactorization, A::AbstractSparseArray{Float64, Int}, b, u, Pl, Pr,
145-
maxiters::Int, abstol,
146-
reltol,
147-
verbose::Bool, assumptions::OperatorAssumptions)
148-
PREALLOCATED_UMFPACK
149-
end
150-
151162
function LinearSolve.init_cacheval(
152163
alg::UMFPACKFactorization, A::LinearSolve.GPUArraysCore.AnyGPUArray, b, u,
153164
Pl, Pr,
@@ -156,6 +167,16 @@ function LinearSolve.init_cacheval(
156167
nothing
157168
end
158169

170+
@static if Base.USE_GPL_LIBS
171+
172+
function LinearSolve.init_cacheval(
173+
alg::UMFPACKFactorization, A::AbstractSparseArray{Float64, Int}, b, u, Pl, Pr,
174+
maxiters::Int, abstol,
175+
reltol,
176+
verbose::Bool, assumptions::OperatorAssumptions)
177+
PREALLOCATED_UMFPACK
178+
end
179+
159180
function LinearSolve.init_cacheval(
160181
alg::UMFPACKFactorization, A::AbstractSparseArray{T, Int64}, b, u,
161182
Pl, Pr,
@@ -211,8 +232,14 @@ function SciMLBase.solve!(
211232
end
212233
end
213234

214-
const PREALLOCATED_KLU = KLU.KLUFactorization(SparseMatrixCSC(0, 0, [1], Int[],
215-
Float64[]))
235+
else
236+
237+
function SciMLBase.solve!(
238+
cache::LinearSolve.LinearCache, alg::UMFPACKFactorization; kwargs...)
239+
error("UMFPACKFactorization requires GPL libraries (UMFPACK). Rebuild Julia with USE_GPL_LIBS=1 or use an alternative algorithm like SparspakFactorization")
240+
end
241+
242+
end # @static if Base.USE_GPL_LIBS
216243

217244
function LinearSolve.init_cacheval(
218245
alg::KLUFactorization, A::AbstractArray, b, u, Pl,
@@ -222,14 +249,6 @@ function LinearSolve.init_cacheval(
222249
nothing
223250
end
224251

225-
function LinearSolve.init_cacheval(
226-
alg::KLUFactorization, A::AbstractSparseArray{Float64, Int64}, b, u, Pl, Pr,
227-
maxiters::Int, abstol,
228-
reltol,
229-
verbose::Bool, assumptions::OperatorAssumptions)
230-
PREALLOCATED_KLU
231-
end
232-
233252
function LinearSolve.init_cacheval(
234253
alg::KLUFactorization, A::LinearSolve.GPUArraysCore.AnyGPUArray, b, u,
235254
Pl, Pr,
@@ -238,6 +257,19 @@ function LinearSolve.init_cacheval(
238257
nothing
239258
end
240259

260+
@static if Base.USE_GPL_LIBS
261+
262+
const PREALLOCATED_KLU = KLU.KLUFactorization(SparseMatrixCSC(0, 0, [1], Int[],
263+
Float64[]))
264+
265+
function LinearSolve.init_cacheval(
266+
alg::KLUFactorization, A::AbstractSparseArray{Float64, Int64}, b, u, Pl, Pr,
267+
maxiters::Int, abstol,
268+
reltol,
269+
verbose::Bool, assumptions::OperatorAssumptions)
270+
PREALLOCATED_KLU
271+
end
272+
241273
function LinearSolve.init_cacheval(
242274
alg::KLUFactorization, A::AbstractSparseArray{Float64, Int32}, b, u, Pl, Pr,
243275
maxiters::Int, abstol,
@@ -247,7 +279,6 @@ function LinearSolve.init_cacheval(
247279
0, 0, [Int32(1)], Int32[], Float64[]))
248280
end
249281

250-
# TODO: guard this against errors
251282
function SciMLBase.solve!(cache::LinearSolve.LinearCache, alg::KLUFactorization; kwargs...)
252283
A = cache.A
253284
A = convert(AbstractMatrix, A)
@@ -282,6 +313,24 @@ function SciMLBase.solve!(cache::LinearSolve.LinearCache, alg::KLUFactorization;
282313
end
283314
end
284315

316+
else
317+
318+
function SciMLBase.solve!(cache::LinearSolve.LinearCache, alg::KLUFactorization; kwargs...)
319+
error("KLUFactorization requires GPL libraries (KLU/SuiteSparse). Rebuild Julia with USE_GPL_LIBS=1 or use an alternative algorithm like SparspakFactorization")
320+
end
321+
322+
end # @static if Base.USE_GPL_LIBS
323+
324+
function LinearSolve.init_cacheval(alg::CHOLMODFactorization,
325+
A::AbstractArray, b, u,
326+
Pl, Pr,
327+
maxiters::Int, abstol, reltol,
328+
verbose::Bool, assumptions::OperatorAssumptions)
329+
nothing
330+
end
331+
332+
@static if Base.USE_GPL_LIBS
333+
285334
const PREALLOCATED_CHOLMOD = cholesky(sparse(reshape([1.0], 1, 1)))
286335

287336
function LinearSolve.init_cacheval(alg::CHOLMODFactorization,
@@ -302,13 +351,7 @@ function LinearSolve.init_cacheval(alg::CHOLMODFactorization,
302351
cholesky(sparse(reshape([one(T)], 1, 1)))
303352
end
304353

305-
function LinearSolve.init_cacheval(alg::CHOLMODFactorization,
306-
A::AbstractArray, b, u,
307-
Pl, Pr,
308-
maxiters::Int, abstol, reltol,
309-
verbose::Bool, assumptions::OperatorAssumptions)
310-
nothing
311-
end
354+
end # @static if Base.USE_GPL_LIBS
312355

313356
function LinearSolve.init_cacheval(alg::NormalCholeskyFactorization,
314357
A::Union{AbstractSparseArray{T}, LinearSolve.GPUArraysCore.AnyGPUArray,
@@ -321,39 +364,58 @@ end
321364
# Specialize QR for the non-square case
322365
# Missing ldiv! definitions: https://github.com/JuliaSparse/SparseArrays.jl/issues/242
323366
function LinearSolve._ldiv!(x::Vector,
324-
A::Union{QR, LinearAlgebra.QRCompactWY,
325-
SparseArrays.SPQR.QRSparse,
326-
SparseArrays.CHOLMOD.Factor}, b::Vector)
367+
A::Union{QR, LinearAlgebra.QRCompactWY}, b::Vector)
327368
x .= A \ b
328369
end
329370

330371
function LinearSolve._ldiv!(x::AbstractVector,
331-
A::Union{QR, LinearAlgebra.QRCompactWY,
332-
SparseArrays.SPQR.QRSparse,
333-
SparseArrays.CHOLMOD.Factor}, b::AbstractVector)
372+
A::Union{QR, LinearAlgebra.QRCompactWY}, b::AbstractVector)
334373
x .= A \ b
335374
end
336375

337376
# Ambiguity removal
338377
function LinearSolve._ldiv!(::SVector,
339-
A::Union{SparseArrays.CHOLMOD.Factor, LinearAlgebra.QR,
340-
LinearAlgebra.QRCompactWY, SparseArrays.SPQR.QRSparse},
378+
A::Union{LinearAlgebra.QR, LinearAlgebra.QRCompactWY},
379+
b::AbstractVector)
380+
(A \ b)
381+
end
382+
function LinearSolve._ldiv!(::SVector,
383+
A::Union{LinearAlgebra.QR, LinearAlgebra.QRCompactWY},
384+
b::SVector)
385+
(A \ b)
386+
end
387+
388+
@static if Base.USE_GPL_LIBS
389+
# SPQR and CHOLMOD Factor support
390+
function LinearSolve._ldiv!(x::Vector,
391+
A::Union{SparseArrays.SPQR.QRSparse, SparseArrays.CHOLMOD.Factor}, b::Vector)
392+
x .= A \ b
393+
end
394+
395+
function LinearSolve._ldiv!(x::AbstractVector,
396+
A::Union{SparseArrays.SPQR.QRSparse, SparseArrays.CHOLMOD.Factor}, b::AbstractVector)
397+
x .= A \ b
398+
end
399+
400+
function LinearSolve._ldiv!(::SVector,
401+
A::Union{SparseArrays.CHOLMOD.Factor, SparseArrays.SPQR.QRSparse},
341402
b::AbstractVector)
342403
(A \ b)
343404
end
344405
function LinearSolve._ldiv!(::SVector,
345-
A::Union{SparseArrays.CHOLMOD.Factor, LinearAlgebra.QR,
346-
LinearAlgebra.QRCompactWY, SparseArrays.SPQR.QRSparse},
406+
A::Union{SparseArrays.CHOLMOD.Factor, SparseArrays.SPQR.QRSparse},
347407
b::SVector)
348408
(A \ b)
349409
end
410+
end # @static if Base.USE_GPL_LIBS
350411

351412
function LinearSolve.pattern_changed(fact, A::SparseArrays.SparseMatrixCSC)
352413
!(SparseArrays.decrement(SparseArrays.getcolptr(A)) ==
353414
fact.colptr && SparseArrays.decrement(SparseArrays.getrowval(A)) ==
354415
fact.rowval)
355416
end
356417

418+
@static if Base.USE_GPL_LIBS
357419
function LinearSolve.defaultalg(
358420
A::AbstractSparseMatrixCSC{<:Union{Float64, ComplexF64}, Ti}, b,
359421
assump::OperatorAssumptions{Bool}) where {Ti}
@@ -367,6 +429,22 @@ function LinearSolve.defaultalg(
367429
LinearSolve.DefaultLinearSolver(LinearSolve.DefaultAlgorithmChoice.QRFactorization)
368430
end
369431
end
432+
else
433+
function LinearSolve.defaultalg(
434+
A::AbstractSparseMatrixCSC{<:Union{Float64, ComplexF64}, Ti}, b,
435+
assump::OperatorAssumptions{Bool}) where {Ti}
436+
ext = Base.get_extension(LinearSolve, :LinearSolveSparspakExt)
437+
if assump.issq && ext !== nothing
438+
LinearSolve.DefaultLinearSolver(LinearSolve.DefaultAlgorithmChoice.SparspakFactorization)
439+
elseif !assump.issq
440+
LinearSolve.DefaultLinearSolver(LinearSolve.DefaultAlgorithmChoice.QRFactorization)
441+
elseif ext === nothing
442+
error("SparspakFactorization required for sparse matrix LU without GPL libraries. Do `using Sparspak` to enable this functionality")
443+
else
444+
error("Unreachable reached. Please report this error with a reproducer.")
445+
end
446+
end
447+
end # @static if Base.USE_GPL_LIBS
370448

371449
# SPQR Handling
372450
function LinearSolve.init_cacheval(

src/default.jl

Lines changed: 15 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -396,9 +396,17 @@ function algchoice_to_alg(alg::Symbol)
396396
elseif alg === :SparspakFactorization
397397
SparspakFactorization(throwerror = false)
398398
elseif alg === :KLUFactorization
399-
KLUFactorization()
399+
@static if Base.USE_GPL_LIBS
400+
KLUFactorization()
401+
else
402+
error("KLUFactorization requires GPL libraries. Rebuild Julia with USE_GPL_LIBS=1 or use a different algorithm")
403+
end
400404
elseif alg === :UMFPACKFactorization
401-
UMFPACKFactorization()
405+
@static if Base.USE_GPL_LIBS
406+
UMFPACKFactorization()
407+
else
408+
error("UMFPACKFactorization requires GPL libraries. Rebuild Julia with USE_GPL_LIBS=1 or use a different algorithm")
409+
end
402410
elseif alg === :KrylovJL_GMRES
403411
KrylovJL_GMRES()
404412
elseif alg === :GenericLUFactorization
@@ -408,7 +416,11 @@ function algchoice_to_alg(alg::Symbol)
408416
elseif alg === :BunchKaufmanFactorization
409417
BunchKaufmanFactorization()
410418
elseif alg === :CHOLMODFactorization
411-
CHOLMODFactorization()
419+
@static if Base.USE_GPL_LIBS
420+
CHOLMODFactorization()
421+
else
422+
error("CHOLMODFactorization requires GPL libraries. Rebuild Julia with USE_GPL_LIBS=1 or use CholeskyFactorization instead")
423+
end
412424
elseif alg === :CholeskyFactorization
413425
CholeskyFactorization()
414426
elseif alg === :NormalCholeskyFactorization

0 commit comments

Comments
 (0)