diff --git a/stdlib/LinearAlgebra/src/bidiag.jl b/stdlib/LinearAlgebra/src/bidiag.jl index 1129b1e25cc47..576d9c9443359 100644 --- a/stdlib/LinearAlgebra/src/bidiag.jl +++ b/stdlib/LinearAlgebra/src/bidiag.jl @@ -202,6 +202,18 @@ function svd(M::Bidiagonal; full::Bool = false) svd!(copy(M), full = full) end +################### +# opnorm # +################### + +function _opnorm1Inf(B::Bidiagonal, p) + size(B, 1) == 1 && return norm(first(B.dv)) + case = xor(p == 1, istriu(B)) + normd1, normdend = norm(first(B.dv)), norm(last(B.dv)) + normd1, normdend = case ? (zero(normd1), normdend) : (normd1, zero(normdend)) + return max(mapreduce(t -> sum(norm, t), max, zip(view(B.dv, (1:length(B.ev)) .+ !case), B.ev)), normdend) +end + #################### # Generic routines # #################### diff --git a/stdlib/LinearAlgebra/src/diagonal.jl b/stdlib/LinearAlgebra/src/diagonal.jl index 69c74e23b07d6..77c6e39330b86 100644 --- a/stdlib/LinearAlgebra/src/diagonal.jl +++ b/stdlib/LinearAlgebra/src/diagonal.jl @@ -548,6 +548,16 @@ function svd(D::Diagonal{<:Number}) return SVD(Up, S[piv], copy(Vp')) end +_opnorm12Inf(A::Diagonal, p) = maximum(opnorm(x, p) for x in A.diag) + +function opnorm(A::Diagonal, p::Real=2) + if p == 1 || p == Inf || p == 2 + return _opnorm12Inf(A, p) + else + throw(ArgumentError("invalid p-norm p=$p. Valid: 1, 2, Inf")) + end +end + # disambiguation methods: * of Diagonal and Adj/Trans AbsVec *(x::Adjoint{<:Any,<:AbstractVector}, D::Diagonal) = Adjoint(map((t,s) -> t'*s, D.diag, parent(x))) *(x::Adjoint{<:Any,<:AbstractVector}, D::Diagonal, y::AbstractVector) = diff --git a/stdlib/LinearAlgebra/src/special.jl b/stdlib/LinearAlgebra/src/special.jl index 4043762f75afa..e20d7147a43a0 100644 --- a/stdlib/LinearAlgebra/src/special.jl +++ b/stdlib/LinearAlgebra/src/special.jl @@ -338,3 +338,16 @@ end ==(A::Bidiagonal, B::SymTridiagonal) = iszero(B.ev) && iszero(A.ev) && A.dv == B.dv ==(B::SymTridiagonal, A::Bidiagonal) = A == B + +# op norm for structured matrices other than diagonal +# the _opnorm1Inf methods are housed in the structured type's respective files + +function opnorm(A::Union{Bidiagonal,Tridiagonal,SymTridiagonal}, p::Real=2) + if p == 2 + return opnorm2(A) + elseif p == 1 || p == Inf + return _opnorm1Inf(A, p) + else + throw(ArgumentError("invalid p-norm p=$p. Valid: 1, 2, Inf")) + end +end diff --git a/stdlib/LinearAlgebra/src/tridiag.jl b/stdlib/LinearAlgebra/src/tridiag.jl index 46b0d02a46711..eb01e1d8efbd8 100644 --- a/stdlib/LinearAlgebra/src/tridiag.jl +++ b/stdlib/LinearAlgebra/src/tridiag.jl @@ -647,5 +647,40 @@ function SymTridiagonal{T}(M::Tridiagonal) where T end end +################### +# opnorms # +################### + +# Tridiagonal + +function _opnorm1Inf(A::Tridiagonal, p) + size(A, 1) == 1 && return norm(first(A.d)) + case = p == Inf + lowerrange, upperrange = case ? (1:length(A.dl)-1, 2:length(A.dl)) : (2:length(A.dl), 1:length(A.dl)-1) + normfirst, normend = case ? (norm(first(A.d))+norm(first(A.du)), norm(last(A.dl))+norm(last(A.d))) : (norm(first(A.d))+norm(first(A.dl)), norm(last(A.du))+norm(last(A.d))) + + return max( + mapreduce(t -> sum(norm, t), + max, + zip(view(A.d, (2:length(A.d)-1)), view(A.dl, lowerrange), view(A.du, upperrange)) + ), + normfirst, normend) +end + +# SymTridiagonal + +function _opnorm1Inf(A::SymTridiagonal, p::Real) + size(A, 1) == 1 && return norm(first(A.dv)) + lowerrange, upperrange = 1:length(A.ev)-1, 2:length(A.ev) + normfirst, normend = norm(first(A.dv))+norm(first(A.ev)), norm(last(A.ev))+norm(last(A.dv)) + + return max( + mapreduce(t -> sum(norm, t), + max, + zip(view(A.dv, (2:length(A.dv)-1)), view(A.ev, lowerrange), view(A.ev, upperrange)) + ), + normfirst, normend) +end + Base._sum(A::Tridiagonal, ::Colon) = sum(A.d) + sum(A.dl) + sum(A.du) Base._sum(A::SymTridiagonal, ::Colon) = sum(A.dv) + 2sum(A.ev) diff --git a/stdlib/LinearAlgebra/test/bidiag.jl b/stdlib/LinearAlgebra/test/bidiag.jl index 250bbbc6cb2ee..ba4789da79de3 100644 --- a/stdlib/LinearAlgebra/test/bidiag.jl +++ b/stdlib/LinearAlgebra/test/bidiag.jl @@ -409,6 +409,32 @@ end @test vcat((Aub\bb)...) ≈ UpperTriangular(A)\b end +@testset "opnorms" begin + B = Bidiagonal([1,-2,3,-4], [1,2,3], 'U') + + @test opnorm(B, 1) == opnorm(Matrix(B), 1) + @test opnorm(B, 2) ≈ opnorm(Matrix(B), 2) + @test opnorm(B, Inf) == opnorm(Matrix(B), Inf) + + B = Bidiagonal([1,-2,3,-4], [1,2,3], 'L') + + @test opnorm(B, 1) == opnorm(Matrix(B), 1) + @test opnorm(B, 2) ≈ opnorm(Matrix(B), 2) + @test opnorm(B, Inf) == opnorm(Matrix(B), Inf) + + B = Bidiagonal([2], Int[], 'L') + + @test opnorm(B, 1) == opnorm(Matrix(B), 1) + @test opnorm(B, 2) ≈ opnorm(Matrix(B), 2) + @test opnorm(B, Inf) == opnorm(Matrix(B), Inf) + + B = Bidiagonal([2], Int[], 'U') + + @test opnorm(B, 1) == opnorm(Matrix(B), 1) + @test opnorm(B, 2) ≈ opnorm(Matrix(B), 2) + @test opnorm(B, Inf) == opnorm(Matrix(B), Inf) +end + @testset "sum" begin @test sum(Bidiagonal([1,2,3], [1,2], :U)) == 9 @test sum(Bidiagonal([1,2,3], [1,2], :L)) == 9 diff --git a/stdlib/LinearAlgebra/test/diagonal.jl b/stdlib/LinearAlgebra/test/diagonal.jl index 69bebeffb2ff2..a5a97a3780520 100644 --- a/stdlib/LinearAlgebra/test/diagonal.jl +++ b/stdlib/LinearAlgebra/test/diagonal.jl @@ -539,6 +539,26 @@ end end end +@testset "opnorms" begin + D = Diagonal([1,-2,3,-4]) + + @test opnorm(D, 1) == opnorm(Matrix(D), 1) + @test opnorm(D, 2) ≈ opnorm(Matrix(D), 2) + @test opnorm(D, Inf) == opnorm(Matrix(D), Inf) + + D = Diagonal([-11]) + @test opnorm(D, 1) == opnorm(Matrix(D), 1) + @test opnorm(D, 2) ≈ opnorm(Matrix(D), 2) + @test opnorm(D, Inf) == opnorm(Matrix(D), Inf) + + # block diagonal matrices + D = Diagonal([[1 2; 3 4], [5 6; 7 8]]) + A = [1 2 0 0; 3 4 0 0; 0 0 5 6; 0 0 7 8] # full matrix of D + @test opnorm(D, 1) == opnorm(A, 1) + @test opnorm(D, 2) ≈ opnorm(A, 2) + @test opnorm(D, Inf) == opnorm(A, Inf) +end + @testset "eigenvalue sorting" begin D = Diagonal([0.4, 0.2, -1.3]) @test eigvals(D) == eigen(D).values == [0.4, 0.2, -1.3] # not sorted by default diff --git a/stdlib/LinearAlgebra/test/tridiag.jl b/stdlib/LinearAlgebra/test/tridiag.jl index 2706b4c768b65..13e05892c329d 100644 --- a/stdlib/LinearAlgebra/test/tridiag.jl +++ b/stdlib/LinearAlgebra/test/tridiag.jl @@ -437,6 +437,30 @@ end @test A90\b90 ≈ inv(A90)*b90 end +@testset "opnorms" begin + T = Tridiagonal([1,2,3], [1,-2,3,-4], [1,2,3]) + + @test opnorm(T, 1) == opnorm(Matrix(T), 1) + @test_skip opnorm(T, 2) ≈ opnorm(Matrix(T), 2) # currently missing + @test opnorm(T, Inf) == opnorm(Matrix(T), Inf) + + S = SymTridiagonal([1,-2,3,-4], [1,2,3]) + + @test opnorm(S, 1) == opnorm(Matrix(S), 1) + @test_skip opnorm(S, 2) ≈ opnorm(Matrix(S), 2) # currently missing + @test opnorm(S, Inf) == opnorm(Matrix(S), Inf) + + T = Tridiagonal(Int[], [-5], Int[]) + @test opnorm(T, 1) == opnorm(Matrix(T), 1) + @test_skip opnorm(T, 2) ≈ opnorm(Matrix(T), 2) # currently missing + @test opnorm(T, Inf) == opnorm(Matrix(T), Inf) + + S = SymTridiagonal(T) + @test opnorm(S, 1) == opnorm(Matrix(S), 1) + @test_skip opnorm(S, 2) ≈ opnorm(Matrix(S), 2) # currently missing + @test opnorm(S, Inf) == opnorm(Matrix(S), Inf) +end + @testset "singular values of SymTridiag" begin @test svdvals(SymTridiagonal([-4,2,3], [0,0])) ≈ [4,3,2] @test svdvals(SymTridiagonal(collect(0.:10.), zeros(10))) ≈ reverse(0:10)