-
Notifications
You must be signed in to change notification settings - Fork 5
add support for BigFloats via a new extension #87
New issue
Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.
By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.
Already on GitHub? Sign in to your account
Merged
Jutho
merged 26 commits into
QuantumKitHub:main
from
sanderdemeyer:GenericLinearAlgebraExt
Nov 5, 2025
Merged
Changes from 23 commits
Commits
Show all changes
26 commits
Select commit
Hold shift + click to select a range
93f37f5
add GenericExt
sanderdemeyer 5e08f7a
include weak dependencies
sanderdemeyer 92534bc
Update ext/MatrixAlgebraKitGenericExt.jl
sanderdemeyer ddaf055
Update ext/MatrixAlgebraKitGenericExt.jl
sanderdemeyer 88b9902
comments from Lukas
sanderdemeyer 3eee758
remove `BigFloat_LQ_Householder`
sanderdemeyer f8d331e
Change struct names and relax type restrictions
sanderdemeyer 03caa56
Split extensions
sanderdemeyer faa0921
fix copies
sanderdemeyer 9ed02bd
fix type instability
sanderdemeyer c3e069e
Name change
sanderdemeyer 7af4d52
Update ext/MatrixAlgebraKitGenericLinearAlgebraExt.jl
sanderdemeyer c6ba4e4
Update ext/MatrixAlgebraKitGenericLinearAlgebraExt.jl
sanderdemeyer fdb1222
Update ext/MatrixAlgebraKitGenericLinearAlgebraExt.jl
sanderdemeyer 49404bf
Merge branch 'QuantumKitHub:main' into GenericLinearAlgebraExt
sanderdemeyer 837a7a4
Resolve comments
sanderdemeyer 78f92f7
change folder names
sanderdemeyer fb1994f
add docs and namechange
sanderdemeyer 3a36446
Remove unnecessary type parameters
lkdvos e6f64b2
simplify SVD and remove allocations
lkdvos 4323b17
simplify Eigh and remove allocations
lkdvos f99112b
simplify QR
lkdvos 89fd42b
simplify Eig and remove allocations
lkdvos 35af44d
switch to `lmul!`
lkdvos 4fe0fe1
docstring improvements
lkdvos c08dad4
move Diagonal tests
lkdvos File filter
Filter by extension
Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
There are no files selected for viewing
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,108 @@ | ||
| module MatrixAlgebraKitGenericLinearAlgebraExt | ||
|
|
||
| using MatrixAlgebraKit | ||
| using MatrixAlgebraKit: sign_safe, check_input, diagview | ||
| using GenericLinearAlgebra: svd!, svdvals!, eigen!, eigvals!, Hermitian, qr! | ||
| using LinearAlgebra: I, Diagonal, rmul!, lmul!, transpose!, dot | ||
|
|
||
| function MatrixAlgebraKit.default_svd_algorithm(::Type{T}; kwargs...) where {T <: StridedMatrix{<:Union{BigFloat, Complex{BigFloat}}}} | ||
| return GLA_QRIteration() | ||
| end | ||
|
|
||
| for f! in (:svd_compact!, :svd_full!, :svd_vals!) | ||
| @eval MatrixAlgebraKit.initialize_output(::typeof($f!), A::AbstractMatrix, ::GLA_QRIteration) = nothing | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.svd_compact!(A::AbstractMatrix, USVᴴ, ::GLA_QRIteration) | ||
| F = svd!(A) | ||
| U, S, Vᴴ = F.U, Diagonal(F.S), F.Vt | ||
| return MatrixAlgebraKit.gaugefix!(svd_compact!, U, S, Vᴴ, size(A)...) | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.svd_full!(A::AbstractMatrix, USVᴴ, ::GLA_QRIteration) | ||
| F = svd!(A; full = true) | ||
| U, Vᴴ = F.U, F.Vt | ||
| S = MatrixAlgebraKit.zero!(similar(F.S, (size(U, 2), size(Vᴴ, 1)))) | ||
| diagview(S) .= F.S | ||
| return MatrixAlgebraKit.gaugefix!(svd_full!, U, S, Vᴴ, size(A)...) | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.svd_vals!(A::AbstractMatrix, S, ::GLA_QRIteration) | ||
| return svdvals!(A) | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.default_eigh_algorithm(::Type{T}; kwargs...) where {T <: StridedMatrix{<:Union{BigFloat, Complex{BigFloat}}}} | ||
| return GLA_QRIteration(; kwargs...) | ||
| end | ||
|
|
||
| for f! in (:eigh_full!, :eigh_vals!) | ||
| @eval MatrixAlgebraKit.initialize_output(::typeof($f!), A::AbstractMatrix, ::GLA_QRIteration) = nothing | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.eigh_full!(A::AbstractMatrix, DV, ::GLA_QRIteration) | ||
| eigval, eigvec = eigen!(Hermitian(A); sortby = real) | ||
| return Diagonal(eigval::AbstractVector{real(eltype(A))}), eigvec::AbstractMatrix{eltype(A)} | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.eigh_vals!(A::AbstractMatrix, D, ::GLA_QRIteration) | ||
| return eigvals!(Hermitian(A); sortby = real) | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.default_qr_algorithm(::Type{T}; kwargs...) where {T <: StridedMatrix{<:Union{BigFloat, Complex{BigFloat}}}} | ||
| return GLA_HouseholderQR(; kwargs...) | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.qr_full!(A::AbstractMatrix, QR, alg::GLA_HouseholderQR) | ||
| check_input(qr_full!, A, QR, alg) | ||
| Q, R = QR | ||
| return _gla_householder_qr!(A, Q, R; alg.kwargs...) | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.qr_compact!(A::AbstractMatrix, QR, alg::GLA_HouseholderQR) | ||
| check_input(qr_compact!, A, QR, alg) | ||
| Q, R = QR | ||
| return _gla_householder_qr!(A, Q, R; alg.kwargs...) | ||
| end | ||
|
|
||
| function _gla_householder_qr!(A::AbstractMatrix, Q, R; positive = false, blocksize = 1, pivoted = false) | ||
| pivoted && throw(ArgumentError("Only pivoted = false implemented for GLA_HouseholderQR.")) | ||
| (blocksize == 1) || throw(ArgumentError("Only blocksize = 1 implemented for GLA_HouseholderQR.")) | ||
|
|
||
| m, n = size(A) | ||
| k = min(m, n) | ||
| computeR = length(R) > 0 | ||
| compact = k < m | ||
| Q̃, R̃ = qr!(A) | ||
|
|
||
| Q = compact ? copyto!(Q, Q̃) : rmul!(MatrixAlgebraKit.one!(Q), Q̃) | ||
| if positive | ||
| @inbounds for j in 1:k | ||
| s = sign_safe(R̃[j, j]) | ||
| @simd for i in 1:m | ||
| Q[i, j] *= s | ||
| end | ||
| end | ||
| end | ||
| if computeR | ||
| if positive | ||
| @inbounds for j in n:-1:1 | ||
| @simd for i in 1:min(k, j) | ||
| R[i, j] = R̃[i, j] * conj(sign_safe(R̃[i, i])) | ||
| end | ||
| @simd for i in (min(k, j) + 1):size(R, 1) | ||
| R[i, j] = zero(eltype(R)) | ||
| end | ||
| end | ||
| else | ||
| R[1:k, :] .= R̃ | ||
| MatrixAlgebraKit.zero!(@view(R[(k + 1):end, :])) | ||
| end | ||
| end | ||
| return Q, R | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.default_lq_algorithm(::Type{T}; kwargs...) where {T <: StridedMatrix{<:Union{BigFloat, Complex{BigFloat}}}} | ||
| return MatrixAlgebraKit.LQViaTransposedQR(GLA_HouseholderQR(; kwargs...)) | ||
| end | ||
|
|
||
| end | ||
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,25 @@ | ||
| module MatrixAlgebraKitGenericSchurExt | ||
|
|
||
| using MatrixAlgebraKit | ||
| using MatrixAlgebraKit: check_input | ||
| using LinearAlgebra: Diagonal | ||
| using GenericSchur | ||
|
|
||
| function MatrixAlgebraKit.default_eig_algorithm(::Type{T}; kwargs...) where {T <: StridedMatrix{<:Union{BigFloat, Complex{BigFloat}}}} | ||
| return GS_QRIteration(; kwargs...) | ||
| end | ||
|
|
||
| for f! in (:eig_full!, :eig_vals!) | ||
| @eval MatrixAlgebraKit.initialize_output(::typeof($f!), A::AbstractMatrix, ::GS_QRIteration) = nothing | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.eig_full!(A::AbstractMatrix, DV, ::GS_QRIteration) | ||
| D, V = GenericSchur.eigen!(A) | ||
| return Diagonal(D), V | ||
| end | ||
|
|
||
| function MatrixAlgebraKit.eig_vals!(A::AbstractMatrix, D, ::GS_QRIteration) | ||
| return GenericSchur.eigvals!(A) | ||
| end | ||
|
|
||
| end |
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
This file contains hidden or bidirectional Unicode text that may be interpreted or compiled differently than what appears below. To review, open the file in an editor that reveals hidden Unicode characters.
Learn more about bidirectional Unicode characters
| Original file line number | Diff line number | Diff line change |
|---|---|---|
| @@ -0,0 +1,118 @@ | ||
| using MatrixAlgebraKit | ||
| using Test | ||
| using TestExtras | ||
| using StableRNGs | ||
| using LinearAlgebra: LinearAlgebra, Diagonal, I | ||
| using MatrixAlgebraKit: TruncatedAlgorithm, diagview, norm | ||
| using GenericLinearAlgebra | ||
|
|
||
| const eltypes = (BigFloat, Complex{BigFloat}) | ||
|
|
||
| @testset "eigh_full! for T = $T" for T in eltypes | ||
| rng = StableRNG(123) | ||
| m = 54 | ||
| alg = GLA_QRIteration() | ||
|
|
||
| A = randn(rng, T, m, m) | ||
| A = (A + A') / 2 | ||
|
|
||
| D, V = @constinferred eigh_full(A; alg) | ||
| @test A * V ≈ V * D | ||
| @test isunitary(V) | ||
| @test all(isreal, D) | ||
|
|
||
| D2, V2 = eigh_full!(copy(A), (D, V), alg) | ||
| @test D2 ≈ D | ||
| @test V2 ≈ V | ||
|
|
||
| D3 = @constinferred eigh_vals(A, alg) | ||
| @test D ≈ Diagonal(D3) | ||
| end | ||
|
|
||
| @testset "eigh_trunc! for T = $T" for T in eltypes | ||
| rng = StableRNG(123) | ||
| m = 54 | ||
| alg = GLA_QRIteration() | ||
| A = randn(rng, T, m, m) | ||
| A = A * A' | ||
| A = (A + A') / 2 | ||
| Ac = similar(A) | ||
| D₀ = reverse(eigh_vals(A)) | ||
|
|
||
| r = m - 2 | ||
| s = 1 + sqrt(eps(real(T))) | ||
| atol = sqrt(eps(real(T))) | ||
|
|
||
| D1, V1, ϵ1 = @constinferred eigh_trunc(A; alg, trunc = truncrank(r)) | ||
| Dfull, Vfull = eigh_full(A; alg) | ||
| @test length(diagview(D1)) == r | ||
| @test isisometric(V1) | ||
| @test A * V1 ≈ V1 * D1 | ||
| @test LinearAlgebra.opnorm(A - V1 * D1 * V1') ≈ D₀[r + 1] | ||
| @test ϵ1 ≈ norm(view(D₀, (r + 1):m)) atol = atol | ||
|
|
||
| trunc = trunctol(; atol = s * D₀[r + 1]) | ||
| D2, V2, ϵ2 = @constinferred eigh_trunc(A; alg, trunc) | ||
| @test length(diagview(D2)) == r | ||
| @test isisometric(V2) | ||
| @test A * V2 ≈ V2 * D2 | ||
| @test ϵ2 ≈ norm(view(D₀, (r + 1):m)) atol = atol | ||
|
|
||
| s = 1 - sqrt(eps(real(T))) | ||
| trunc = truncerror(; atol = s * norm(@view(D₀[r:end]), 1), p = 1) | ||
| D3, V3, ϵ3 = @constinferred eigh_trunc(A; alg, trunc) | ||
| @test length(diagview(D3)) == r | ||
| @test A * V3 ≈ V3 * D3 | ||
| @test ϵ3 ≈ norm(view(D₀, (r + 1):m)) atol = atol | ||
|
|
||
| # test for same subspace | ||
| @test V1 * (V1' * V2) ≈ V2 | ||
| @test V2 * (V2' * V1) ≈ V1 | ||
| @test V1 * (V1' * V3) ≈ V3 | ||
| @test V3 * (V3' * V1) ≈ V1 | ||
| end | ||
|
|
||
| @testset "eigh_trunc! specify truncation algorithm T = $T" for T in eltypes | ||
| rng = StableRNG(123) | ||
| m = 4 | ||
| atol = sqrt(eps(real(T))) | ||
| V = qr_compact(randn(rng, T, m, m))[1] | ||
| D = Diagonal(real(T)[0.9, 0.3, 0.1, 0.01]) | ||
| A = V * D * V' | ||
| A = (A + A') / 2 | ||
| alg = TruncatedAlgorithm(GLA_QRIteration(), truncrank(2)) | ||
| D2, V2, ϵ2 = @constinferred eigh_trunc(A; alg) | ||
| @test diagview(D2) ≈ diagview(D)[1:2] | ||
| @test_throws ArgumentError eigh_trunc(A; alg, trunc = (; maxrank = 2)) | ||
| @test ϵ2 ≈ norm(diagview(D)[3:4]) atol = atol | ||
|
|
||
| alg = TruncatedAlgorithm(GLA_QRIteration(), truncerror(; atol = 0.2)) | ||
| D3, V3, ϵ3 = @constinferred eigh_trunc(A; alg) | ||
| @test diagview(D3) ≈ diagview(D)[1:2] | ||
| @test ϵ3 ≈ norm(diagview(D)[3:4]) atol = atol | ||
| end | ||
|
|
||
| @testset "eigh for Diagonal{$T}" for T in eltypes | ||
lkdvos marked this conversation as resolved.
Outdated
Show resolved
Hide resolved
|
||
| rng = StableRNG(123) | ||
| m = 54 | ||
| Ad = randn(rng, T, m) | ||
| Ad .+= conj.(Ad) | ||
| A = Diagonal(Ad) | ||
| atol = sqrt(eps(real(T))) | ||
|
|
||
| D, V = @constinferred eigh_full(A) | ||
| @test D isa Diagonal{real(T)} && size(D) == size(A) | ||
| @test V isa Diagonal{T} && size(V) == size(A) | ||
| @test A * V ≈ V * D | ||
|
|
||
| D2 = @constinferred eigh_vals(A) | ||
| @test D2 isa AbstractVector{real(T)} && length(D2) == m | ||
| @test diagview(D) ≈ D2 | ||
|
|
||
| A2 = Diagonal(T[0.9, 0.3, 0.1, 0.01]) | ||
| alg = TruncatedAlgorithm(DiagonalAlgorithm(), truncrank(2)) | ||
| D2, V2, ϵ2 = @constinferred eigh_trunc(A2; alg) | ||
| @test diagview(D2) ≈ diagview(A2)[1:2] | ||
| @test ϵ2 ≈ norm(diagview(A2)[3:4]) atol = atol | ||
|
|
||
| end | ||
Oops, something went wrong.
Add this suggestion to a batch that can be applied as a single commit.
This suggestion is invalid because no changes were made to the code.
Suggestions cannot be applied while the pull request is closed.
Suggestions cannot be applied while viewing a subset of changes.
Only one suggestion per line can be applied in a batch.
Add this suggestion to a batch that can be applied as a single commit.
Applying suggestions on deleted lines is not supported.
You must change the existing code in this line in order to create a valid suggestion.
Outdated suggestions cannot be applied.
This suggestion has been applied or marked resolved.
Suggestions cannot be applied from pending reviews.
Suggestions cannot be applied on multi-line comments.
Suggestions cannot be applied while the pull request is queued to merge.
Suggestion cannot be applied right now. Please check back later.
Uh oh!
There was an error while loading. Please reload this page.