Skip to content

Commit 3cf2634

Browse files
committed
parametrize Tridiagonal on the wrapped vector type
1 parent 1323561 commit 3cf2634

File tree

8 files changed

+90
-83
lines changed

8 files changed

+90
-83
lines changed

NEWS.md

+12-9
Original file line numberDiff line numberDiff line change
@@ -121,10 +121,11 @@ This section lists changes that do not have deprecation warnings.
121121
longer present. Use `first(R)` and `last(R)` to obtain
122122
start/stop. ([#20974])
123123

124-
* The `Diagonal`, `Bidiagonal` and `SymTridiagonal` type definitions have changed from
125-
`Diagonal{T}`, `Bidiagonal{T}` and `SymTridiagonal{T}` to `Diagonal{T,V<:AbstractVector{T}}`,
126-
`Bidiagonal{T,V<:AbstractVector{T}}` and `SymTridiagonal{T,V<:AbstractVector{T}}`
127-
respectively ([#22718], [#22925], [#23035]).
124+
* The `Diagonal`, `Bidiagonal`, `Tridiagonal` and `SymTridiagonal` type definitions have
125+
changed from `Diagonal{T}`, `Bidiagonal{T}`, `Tridiagonal{T}` and `SymTridiagonal{T}`
126+
to `Diagonal{T,V<:AbstractVector{T}}`, `Bidiagonal{T,V<:AbstractVector{T}}`,
127+
`Tridiagonal{T,V<:AbstractVector{T}}` and `SymTridiagonal{T,V<:AbstractVector{T}}`
128+
respectively ([#22718], [#22925], [#23035], [#23154]).
128129

129130
* `isapprox(x,y)` now tests `norm(x-y) <= max(atol, rtol*max(norm(x), norm(y)))`
130131
rather than `norm(x-y) <= atol + ...`, and `rtol` defaults to zero
@@ -196,9 +197,10 @@ Library improvements
196197

197198
* `Char`s can now be concatenated with `String`s and/or other `Char`s using `*` ([#22532]).
198199

199-
* `Diagonal`, `Bidiagonal` and `SymTridiagonal` are now parameterized on the type of the
200-
wrapped vectors, allowing `Diagonal`, `Bidiagonal` and `SymTridiagonal` matrices with
201-
arbitrary `AbstractVector`s ([#22718], [#22925], [#23035]).
200+
* `Diagonal`, `Bidiagonal`, `Tridiagonal` and `SymTridiagonal` are now parameterized on
201+
the type of the wrapped vectors, allowing `Diagonal`, `Bidiagonal`, `Tridiagonal` and
202+
`SymTridiagonal` matrices with arbitrary `AbstractVector`s
203+
([#22718], [#22925], [#23035], [#23154]).
202204

203205
* Mutating versions of `randperm` and `randcycle` have been added:
204206
`randperm!` and `randcycle!` ([#22723]).
@@ -261,8 +263,9 @@ Deprecated or removed
261263
* `Bidiagonal` constructors now use a `Symbol` (`:U` or `:L`) for the upper/lower
262264
argument, instead of a `Bool` or a `Char` ([#22703]).
263265

264-
* `Bidiagonal` and `SymTridiagonal` constructors that automatically converted the input
265-
vectors to the same type are deprecated in favor of explicit conversion ([#22925], [#23035]).
266+
* `Bidiagonal`, `Tridiagonal` and `SymTridiagonal` constructors that automatically
267+
converted the input vectors to the same type are deprecated in favor of explicit
268+
conversion ([#22925], [#23035], [#23154].
266269

267270
* Calling `nfields` on a type to find out how many fields its instances have is deprecated.
268271
Use `fieldcount` instead. Use `nfields` only to get the number of fields in a specific object ([#22350]).

base/deprecated.jl

+8
Original file line numberDiff line numberDiff line change
@@ -1652,6 +1652,14 @@ function SymTridiagonal(dv::AbstractVector{T}, ev::AbstractVector{S}) where {T,S
16521652
SymTridiagonal(convert(Vector{R}, dv), convert(Vector{R}, ev))
16531653
end
16541654

1655+
# PR #23154
1656+
# also uncomment constructor tests in test/linalg/tridiag.jl
1657+
function Tridiagonal(dl::AbstractVector{Tl}, d::AbstractVector{Td}, du::AbstractVector{Tu}) where {Tl,Td,Tu}
1658+
depwarn(string("Tridiagonal(dl::AbstractVector{Tl}, d::AbstractVector{Td}, du::AbstractVector{Tu}) ",
1659+
"where {Tl, Td, Tu} is deprecated; convert all vectors to the same type instead."), :Tridiagonal)
1660+
Tridiagonal(map(v->convert(Vector{promote_type(Tl,Td,Tu)}, v), (dl, d, du))...)
1661+
end
1662+
16551663
# PR #23092
16561664
@eval LibGit2 begin
16571665
function prompt(msg::AbstractString; default::AbstractString="", password::Bool=false)

base/linalg/bidiag.jl

+4-2
Original file line numberDiff line numberDiff line change
@@ -153,8 +153,10 @@ promote_rule(::Type{Matrix{T}}, ::Type{Bidiagonal{S}}) where {T,S} = Matrix{prom
153153
#Converting from Bidiagonal to Tridiagonal
154154
Tridiagonal(M::Bidiagonal{T}) where {T} = convert(Tridiagonal{T}, M)
155155
function convert(::Type{Tridiagonal{T}}, A::Bidiagonal) where T
156-
z = zeros(T, size(A)[1]-1)
157-
A.uplo == 'U' ? Tridiagonal(z, convert(Vector{T},A.dv), convert(Vector{T},A.ev)) : Tridiagonal(convert(Vector{T},A.ev), convert(Vector{T},A.dv), z)
156+
dv = convert(AbstractVector{T}, A.dv)
157+
ev = convert(AbstractVector{T}, A.ev)
158+
z = fill!(similar(ev), zero(T))
159+
A.uplo == 'U' ? Tridiagonal(z, dv, ev) : Tridiagonal(ev, dv, z)
158160
end
159161
promote_rule(::Type{Tridiagonal{T}}, ::Type{Bidiagonal{S}}) where {T,S} = Tridiagonal{promote_type(T,S)}
160162

base/linalg/bunchkaufman.jl

+1-1
Original file line numberDiff line numberDiff line change
@@ -116,7 +116,7 @@ julia> A = [1 2 3; 2 1 2; 3 2 1]
116116
julia> F = bkfact(Symmetric(A, :L))
117117
Base.LinAlg.BunchKaufman{Float64,Array{Float64,2}}
118118
D factor:
119-
3×3 Tridiagonal{Float64}:
119+
3×3 Tridiagonal{Float64,Array{Float64,1}}:
120120
1.0 3.0 ⋅
121121
3.0 1.0 0.0
122122
⋅ 0.0 -1.0

base/linalg/lu.jl

+14-13
Original file line numberDiff line numberDiff line change
@@ -327,14 +327,14 @@ end
327327
# Tridiagonal
328328

329329
# See dgttrf.f
330-
function lufact!(A::Tridiagonal{T}, pivot::Union{Val{false}, Val{true}} = Val(true)) where T
330+
function lufact!(A::Tridiagonal{T,V}, pivot::Union{Val{false}, Val{true}} = Val(true)) where {T,V}
331331
n = size(A, 1)
332332
info = 0
333333
ipiv = Vector{BlasInt}(n)
334334
dl = A.dl
335335
d = A.d
336336
du = A.du
337-
du2 = A.du2
337+
du2 = fill!(similar(d, n-2), 0)::V
338338

339339
@inbounds begin
340340
for i = 1:n
@@ -389,12 +389,13 @@ function lufact!(A::Tridiagonal{T}, pivot::Union{Val{false}, Val{true}} = Val(tr
389389
end
390390
end
391391
end
392-
LU{T,Tridiagonal{T}}(A, ipiv, convert(BlasInt, info))
392+
B = Tridiagonal{T,V}(dl, d, du, du2)
393+
LU{T,Tridiagonal{T,V}}(B, ipiv, convert(BlasInt, info))
393394
end
394395

395396
factorize(A::Tridiagonal) = lufact(A)
396397

397-
function getindex(F::Base.LinAlg.LU{T,Tridiagonal{T}}, d::Symbol) where T
398+
function getindex(F::LU{T,Tridiagonal{T,V}}, d::Symbol) where {T,V}
398399
m, n = size(F)
399400
if d == :L
400401
L = Array(Bidiagonal(ones(T, n), F.factors.dl, d))
@@ -419,7 +420,7 @@ function getindex(F::Base.LinAlg.LU{T,Tridiagonal{T}}, d::Symbol) where T
419420
end
420421

421422
# See dgtts2.f
422-
function A_ldiv_B!(A::LU{T,Tridiagonal{T}}, B::AbstractVecOrMat) where T
423+
function A_ldiv_B!(A::LU{T,Tridiagonal{T,V}}, B::AbstractVecOrMat) where {T,V}
423424
n = size(A,1)
424425
if n != size(B,1)
425426
throw(DimensionMismatch("matrix has dimensions ($n,$n) but right hand side has $(size(B,1)) rows"))
@@ -450,7 +451,7 @@ function A_ldiv_B!(A::LU{T,Tridiagonal{T}}, B::AbstractVecOrMat) where T
450451
return B
451452
end
452453

453-
function At_ldiv_B!(A::LU{T,Tridiagonal{T}}, B::AbstractVecOrMat) where T
454+
function At_ldiv_B!(A::LU{T,Tridiagonal{T,V}}, B::AbstractVecOrMat) where {T,V}
454455
n = size(A,1)
455456
if n != size(B,1)
456457
throw(DimensionMismatch("matrix has dimensions ($n,$n) but right hand side has $(size(B,1)) rows"))
@@ -485,7 +486,7 @@ function At_ldiv_B!(A::LU{T,Tridiagonal{T}}, B::AbstractVecOrMat) where T
485486
end
486487

487488
# Ac_ldiv_B!(A::LU{T,Tridiagonal{T}}, B::AbstractVecOrMat) where {T<:Real} = At_ldiv_B!(A,B)
488-
function Ac_ldiv_B!(A::LU{T,Tridiagonal{T}}, B::AbstractVecOrMat) where T
489+
function Ac_ldiv_B!(A::LU{T,Tridiagonal{T,V}}, B::AbstractVecOrMat) where {T,V}
489490
n = size(A,1)
490491
if n != size(B,1)
491492
throw(DimensionMismatch("matrix has dimensions ($n,$n) but right hand side has $(size(B,1)) rows"))
@@ -528,7 +529,7 @@ convert(::Type{Matrix}, F::LU) = convert(Array, convert(AbstractArray, F))
528529
convert(::Type{Array}, F::LU) = convert(Matrix, F)
529530
full(F::LU) = convert(AbstractArray, F)
530531

531-
function convert(::Type{Tridiagonal}, F::Base.LinAlg.LU{T,Tridiagonal{T}}) where T
532+
function convert(::Type{Tridiagonal}, F::Base.LinAlg.LU{T,Tridiagonal{T,V}}) where {T,V}
532533
n = size(F, 1)
533534

534535
dl = copy(F.factors.dl)
@@ -562,12 +563,12 @@ function convert(::Type{Tridiagonal}, F::Base.LinAlg.LU{T,Tridiagonal{T}}) where
562563
end
563564
return Tridiagonal(dl, d, du)
564565
end
565-
convert(::Type{AbstractMatrix}, F::Base.LinAlg.LU{T,Tridiagonal{T}}) where {T} =
566+
convert(::Type{AbstractMatrix}, F::LU{T,Tridiagonal{T,V}}) where {T,V} =
566567
convert(Tridiagonal, F)
567-
convert(::Type{AbstractArray}, F::Base.LinAlg.LU{T,Tridiagonal{T}}) where {T} =
568+
convert(::Type{AbstractArray}, F::LU{T,Tridiagonal{T,V}}) where {T,V} =
568569
convert(AbstractMatrix, F)
569-
convert(::Type{Matrix}, F::Base.LinAlg.LU{T,Tridiagonal{T}}) where {T} =
570+
convert(::Type{Matrix}, F::LU{T,Tridiagonal{T,V}}) where {T,V} =
570571
convert(Array, convert(AbstractArray, F))
571-
convert(::Type{Array}, F::Base.LinAlg.LU{T,Tridiagonal{T}}) where {T} =
572+
convert(::Type{Array}, F::LU{T,Tridiagonal{T,V}}) where {T,V} =
572573
convert(Matrix, F)
573-
full(F::Base.LinAlg.LU{T,Tridiagonal{T}}) where {T} = convert(AbstractArray, F)
574+
full(F::LU{T,Tridiagonal{T,V}}) where {T,V} = convert(AbstractArray, F)

base/linalg/tridiag.jl

+43-57
Original file line numberDiff line numberDiff line change
@@ -399,71 +399,58 @@ function setindex!(A::SymTridiagonal, x, i::Integer, j::Integer)
399399
end
400400

401401
## Tridiagonal matrices ##
402-
struct Tridiagonal{T} <: AbstractMatrix{T}
403-
dl::Vector{T} # sub-diagonal
404-
d::Vector{T} # diagonal
405-
du::Vector{T} # sup-diagonal
406-
du2::Vector{T} # supsup-diagonal for pivoting
402+
struct Tridiagonal{T,V<:AbstractVector{T}} <: AbstractMatrix{T}
403+
dl::V # sub-diagonal
404+
d::V # diagonal
405+
du::V # sup-diagonal
406+
du2::V # supsup-diagonal for pivoting in LU
407+
function Tridiagonal{T}(dl::V, d::V, du::V) where {T,V<:AbstractVector{T}}
408+
n = length(d)
409+
if (length(dl) != n-1 || length(du) != n-1)
410+
throw(ArgumentError(string("cannot construct Tridiagonal from incompatible ",
411+
"lengths of subdiagonal, diagonal and superdiagonal: ",
412+
"($(length(dl)), $(length(d)), $(length(du)))")))
413+
end
414+
new{T,V}(dl, d, du)
415+
end
416+
# constructor used in lufact!
417+
function Tridiagonal{T,V}(dl::V, d::V, du::V, du2::V) where {T,V<:AbstractVector{T}}
418+
new{T,V}(dl, d, du, du2)
419+
end
407420
end
408421

409422
"""
410-
Tridiagonal(dl, d, du)
423+
Tridiagonal(dl::V, d::V, du::V) where V <: AbstractVector
411424
412425
Construct a tridiagonal matrix from the first subdiagonal, diagonal, and first superdiagonal,
413-
respectively. The result is of type `Tridiagonal` and provides efficient specialized linear
426+
respectively. The result is of type `Tridiagonal` and provides efficient specialized linear
414427
solvers, but may be converted into a regular matrix with
415428
[`convert(Array, _)`](@ref) (or `Array(_)` for short).
416429
The lengths of `dl` and `du` must be one less than the length of `d`.
417430
418431
# Examples
419432
```jldoctest
420-
julia> dl = [1; 2; 3]
421-
3-element Array{Int64,1}:
422-
1
423-
2
424-
3
433+
julia> dl = [1, 2, 3];
425434
426-
julia> du = [4; 5; 6]
427-
3-element Array{Int64,1}:
428-
4
429-
5
430-
6
435+
julia> du = [4, 5, 6];
431436
432-
julia> d = [7; 8; 9; 0]
433-
4-element Array{Int64,1}:
434-
7
435-
8
436-
9
437-
0
437+
julia> d = [7, 8, 9, 0];
438438
439439
julia> Tridiagonal(dl, d, du)
440-
4×4 Tridiagonal{Int64}:
440+
4×4 Tridiagonal{Int64,Array{Int64,1}}:
441441
7 4 ⋅ ⋅
442442
1 8 5 ⋅
443443
⋅ 2 9 6
444444
⋅ ⋅ 3 0
445445
```
446446
"""
447-
# Basic constructor takes in three dense vectors of same type
448-
function Tridiagonal(dl::Vector{T}, d::Vector{T}, du::Vector{T}) where T
449-
n = length(d)
450-
if (length(dl) != n-1 || length(du) != n-1)
451-
throw(ArgumentError("cannot make Tridiagonal from incompatible lengths of subdiagonal, diagonal and superdiagonal: ($(length(dl)), $(length(d)), $(length(du))"))
452-
end
453-
Tridiagonal(dl, d, du, zeros(T,n-2))
454-
end
455-
456-
# Construct from diagonals of any abstract vector, any eltype
457-
function Tridiagonal(dl::AbstractVector{Tl}, d::AbstractVector{Td}, du::AbstractVector{Tu}) where {Tl,Td,Tu}
458-
Tridiagonal(map(v->convert(Vector{promote_type(Tl,Td,Tu)}, v), (dl, d, du))...)
459-
end
447+
Tridiagonal(dl::V, d::V, du::V) where {T,V<:AbstractVector{T}} = Tridiagonal{T}(dl, d, du)
460448

461-
# Provide a constructor Tridiagonal(A) similar to the triangulars, diagonal, symmetric
462449
"""
463450
Tridiagonal(A)
464451
465-
returns a `Tridiagonal` array based on (abstract) matrix `A`, using its first lower diagonal,
466-
main diagonal, and first upper diagonal.
452+
Construct a tridiagonal matrix from the first sub-diagonal,
453+
diagonal and first super-diagonal of the matrix `A`.
467454
468455
# Examples
469456
```jldoctest
@@ -475,16 +462,14 @@ julia> A = [1 2 3 4; 1 2 3 4; 1 2 3 4; 1 2 3 4]
475462
1 2 3 4
476463
477464
julia> Tridiagonal(A)
478-
4×4 Tridiagonal{Int64}:
465+
4×4 Tridiagonal{Int64,Array{Int64,1}}:
479466
1 2 ⋅ ⋅
480467
1 2 3 ⋅
481468
⋅ 2 3 4
482469
⋅ ⋅ 3 4
483470
```
484471
"""
485-
function Tridiagonal(A::AbstractMatrix)
486-
return Tridiagonal(diag(A,-1), diag(A), diag(A,+1))
487-
end
472+
Tridiagonal(A::AbstractMatrix) = Tridiagonal(diag(A,-1), diag(A,0), diag(A,1))
488473

489474
size(M::Tridiagonal) = (length(M.d), length(M.d))
490475
function size(M::Tridiagonal, d::Integer)
@@ -512,31 +497,31 @@ convert(::Type{Matrix}, M::Tridiagonal{T}) where {T} = convert(Matrix{T}, M)
512497
convert(::Type{Array}, M::Tridiagonal) = convert(Matrix, M)
513498
full(M::Tridiagonal) = convert(Array, M)
514499
function similar(M::Tridiagonal, ::Type{T}) where T
515-
Tridiagonal{T}(similar(M.dl, T), similar(M.d, T), similar(M.du, T), similar(M.du2, T))
500+
Tridiagonal{T}(similar(M.dl, T), similar(M.d, T), similar(M.du, T))
516501
end
517502

518503
# Operations on Tridiagonal matrices
519-
copy!(dest::Tridiagonal, src::Tridiagonal) = Tridiagonal(copy!(dest.dl, src.dl), copy!(dest.d, src.d), copy!(dest.du, src.du), copy!(dest.du2, src.du2))
504+
copy!(dest::Tridiagonal, src::Tridiagonal) = (copy!(dest.dl, src.dl); copy!(dest.d, src.d); copy!(dest.du, src.du); dest)
520505

521506
#Elementary operations
522-
broadcast(::typeof(abs), M::Tridiagonal) = Tridiagonal(abs.(M.dl), abs.(M.d), abs.(M.du), abs.(M.du2))
523-
broadcast(::typeof(round), M::Tridiagonal) = Tridiagonal(round.(M.dl), round.(M.d), round.(M.du), round.(M.du2))
524-
broadcast(::typeof(trunc), M::Tridiagonal) = Tridiagonal(trunc.(M.dl), trunc.(M.d), trunc.(M.du), trunc.(M.du2))
525-
broadcast(::typeof(floor), M::Tridiagonal) = Tridiagonal(floor.(M.dl), floor.(M.d), floor.(M.du), floor.(M.du2))
526-
broadcast(::typeof(ceil), M::Tridiagonal) = Tridiagonal(ceil.(M.dl), ceil.(M.d), ceil.(M.du), ceil.(M.du2))
507+
broadcast(::typeof(abs), M::Tridiagonal) = Tridiagonal(abs.(M.dl), abs.(M.d), abs.(M.du))
508+
broadcast(::typeof(round), M::Tridiagonal) = Tridiagonal(round.(M.dl), round.(M.d), round.(M.du))
509+
broadcast(::typeof(trunc), M::Tridiagonal) = Tridiagonal(trunc.(M.dl), trunc.(M.d), trunc.(M.du))
510+
broadcast(::typeof(floor), M::Tridiagonal) = Tridiagonal(floor.(M.dl), floor.(M.d), floor.(M.du))
511+
broadcast(::typeof(ceil), M::Tridiagonal) = Tridiagonal(ceil.(M.dl), ceil.(M.d), ceil.(M.du))
527512
for func in (:conj, :copy, :real, :imag)
528513
@eval function ($func)(M::Tridiagonal)
529-
Tridiagonal(($func)(M.dl), ($func)(M.d), ($func)(M.du), ($func)(M.du2))
514+
Tridiagonal(($func)(M.dl), ($func)(M.d), ($func)(M.du))
530515
end
531516
end
532517
broadcast(::typeof(round), ::Type{T}, M::Tridiagonal) where {T<:Integer} =
533-
Tridiagonal(round.(T, M.dl), round.(T, M.d), round.(T, M.du), round.(T, M.du2))
518+
Tridiagonal(round.(T, M.dl), round.(T, M.d), round.(T, M.du))
534519
broadcast(::typeof(trunc), ::Type{T}, M::Tridiagonal) where {T<:Integer} =
535-
Tridiagonal(trunc.(T, M.dl), trunc.(T, M.d), trunc.(T, M.du), trunc.(T, M.du2))
520+
Tridiagonal(trunc.(T, M.dl), trunc.(T, M.d), trunc.(T, M.du))
536521
broadcast(::typeof(floor), ::Type{T}, M::Tridiagonal) where {T<:Integer} =
537-
Tridiagonal(floor.(T, M.dl), floor.(T, M.d), floor.(T, M.du), floor.(T, M.du2))
522+
Tridiagonal(floor.(T, M.dl), floor.(T, M.d), floor.(T, M.du))
538523
broadcast(::typeof(ceil), ::Type{T}, M::Tridiagonal) where {T<:Integer} =
539-
Tridiagonal(ceil.(T, M.dl), ceil.(T, M.d), ceil.(T, M.du), ceil.(T, M.du2))
524+
Tridiagonal(ceil.(T, M.dl), ceil.(T, M.d), ceil.(T, M.du))
540525

541526
transpose(M::Tridiagonal) = Tridiagonal(M.du, M.d, M.dl)
542527
ctranspose(M::Tridiagonal) = conj(transpose(M))
@@ -646,7 +631,8 @@ end
646631
inv(A::Tridiagonal) = inv_usmani(A.dl, A.d, A.du)
647632
det(A::Tridiagonal) = det_usmani(A.dl, A.d, A.du)
648633

649-
convert(::Type{Tridiagonal{T}},M::Tridiagonal) where {T} = Tridiagonal(convert(Vector{T}, M.dl), convert(Vector{T}, M.d), convert(Vector{T}, M.du), convert(Vector{T}, M.du2))
634+
convert(::Type{Tridiagonal{T}},M::Tridiagonal) where {T} =
635+
Tridiagonal(convert(AbstractVector{T}, M.dl), convert(AbstractVector{T}, M.d), convert(AbstractVector{T}, M.du))
650636
convert(::Type{AbstractMatrix{T}},M::Tridiagonal) where {T} = convert(Tridiagonal{T}, M)
651637
convert(::Type{Tridiagonal{T}}, M::SymTridiagonal{T}) where {T} = Tridiagonal(M)
652638
function convert(::Type{SymTridiagonal{T}}, M::Tridiagonal) where T

test/linalg/tridiag.jl

+7
Original file line numberDiff line numberDiff line change
@@ -41,10 +41,17 @@ B = randn(n,2)
4141
@test ST == Matrix(ST)
4242
@test ST.dv === x
4343
@test ST.ev === y
44+
TT = (Tridiagonal(y, x, y))::Tridiagonal{elty, typeof(x)}
45+
@test TT == Matrix(TT)
46+
@test TT.dl === y
47+
@test TT.d === x
48+
@test TT.du === y
4449
end
4550
# enable when deprecations for 0.7 are dropped
4651
# @test_throws MethodError SymTridiagonal(dv, GenericArray(ev))
4752
# @test_throws MethodError SymTridiagonal(GenericArray(dv), ev)
53+
# @test_throws MethodError Tridiagonal(GenericArray(ev), dv, GenericArray(ev))
54+
# @test_throws MethodError Tridiagonal(ev, GenericArray(dv), ev)
4855
end
4956

5057
@testset "size and Array" begin

test/show.jl

+1-1
Original file line numberDiff line numberDiff line change
@@ -598,7 +598,7 @@ A = reshape(1:16,4,4)
598598
@test replstr(Bidiagonal(A,:U)) == "4×4 Bidiagonal{$(Int),Array{$(Int),1}}:\n 1 5 ⋅ ⋅\n ⋅ 6 10 ⋅\n ⋅ ⋅ 11 15\n ⋅ ⋅ ⋅ 16"
599599
@test replstr(Bidiagonal(A,:L)) == "4×4 Bidiagonal{$(Int),Array{$(Int),1}}:\n 1 ⋅ ⋅ ⋅\n 2 6 ⋅ ⋅\n ⋅ 7 11 ⋅\n ⋅ ⋅ 12 16"
600600
@test replstr(SymTridiagonal(A+A')) == "4×4 SymTridiagonal{$(Int),Array{$(Int),1}}:\n 2 7 ⋅ ⋅\n 7 12 17 ⋅\n ⋅ 17 22 27\n ⋅ ⋅ 27 32"
601-
@test replstr(Tridiagonal(diag(A,-1),diag(A),diag(A,+1))) == "4×4 Tridiagonal{$Int}:\n 1 5 ⋅ ⋅\n 2 6 10 ⋅\n ⋅ 7 11 15\n ⋅ ⋅ 12 16"
601+
@test replstr(Tridiagonal(diag(A,-1),diag(A),diag(A,+1))) == "4×4 Tridiagonal{$(Int),Array{$(Int),1}}:\n 1 5 ⋅ ⋅\n 2 6 10 ⋅\n ⋅ 7 11 15\n ⋅ ⋅ 12 16"
602602
@test replstr(UpperTriangular(copy(A))) == "4×4 UpperTriangular{$Int,Array{$Int,2}}:\n 1 5 9 13\n ⋅ 6 10 14\n ⋅ ⋅ 11 15\n ⋅ ⋅ ⋅ 16"
603603
@test replstr(LowerTriangular(copy(A))) == "4×4 LowerTriangular{$Int,Array{$Int,2}}:\n 1 ⋅ ⋅ ⋅\n 2 6 ⋅ ⋅\n 3 7 11 ⋅\n 4 8 12 16"
604604

0 commit comments

Comments
 (0)