Skip to content

Commit

Permalink
Rework some constructors
Browse files Browse the repository at this point in the history
Aside from enforcing 1-indexing, these allow one to coerce some of the
field types via construction (without requiring that the inputs already
have those types). One potential disadvantage is that they may sometimes
copy data.
  • Loading branch information
timholy committed Jul 6, 2018
1 parent 7e891ee commit da19abe
Show file tree
Hide file tree
Showing 10 changed files with 79 additions and 24 deletions.
7 changes: 7 additions & 0 deletions stdlib/LinearAlgebra/src/diagonal.jl
Original file line number Diff line number Diff line change
Expand Up @@ -4,7 +4,14 @@

struct Diagonal{T,V<:AbstractVector{T}} <: AbstractMatrix{T}
diag::V

function Diagonal{T,V}(diag) where {V<:AbstractVector{T}}
@assert is_one_indexed(diag)
new(diag)
end
end
Diagonal(v::AbstractVector{T}) where T = Diagonal{T,typeof(v)}(v)

"""
Diagonal(A::AbstractMatrix)
Expand Down
5 changes: 4 additions & 1 deletion stdlib/LinearAlgebra/src/hessenberg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -63,7 +63,10 @@ hessenberg(A::StridedMatrix{T}) where T =
struct HessenbergQ{T,S<:AbstractMatrix} <: AbstractMatrix{T}
factors::S
τ::Vector{T}
HessenbergQ{T,S}(factors::AbstractMatrix{T}, τ::Vector{T}) where {T,S<:AbstractMatrix} = new(factors, τ)
function HessenbergQ{T,S}(factors, τ) where {T,S<:AbstractMatrix}
@assert is_one_indexed(factors)
new(factors, τ)
end
end
HessenbergQ(factors::AbstractMatrix{T}, τ::Vector{T}) where {T} = HessenbergQ{T,typeof(factors)}(factors, τ)
HessenbergQ(A::Hessenberg) = HessenbergQ(A.factors, A.τ)
Expand Down
8 changes: 7 additions & 1 deletion stdlib/LinearAlgebra/src/ldlt.jl
Original file line number Diff line number Diff line change
@@ -1,8 +1,14 @@
# This file is a part of Julia. License is MIT: https://julialang.org/license

struct LDLt{T,S<:AbstractMatrix} <: Factorization{T}
struct LDLt{T,S<:AbstractMatrix{T}} <: Factorization{T}
data::S

function LDLt{T,S}(data) where {T,S<:AbstractMatrix}
@assert is_one_indexed(data)
new(data)
end
end
LDLt(data::AbstractMatrix{T}) where T = LDLt{T,typeof(data)}(data)

size(S::LDLt) = size(S.data)
size(S::LDLt, i::Integer) = size(S.data, i)
Expand Down
7 changes: 5 additions & 2 deletions stdlib/LinearAlgebra/src/lq.jl
Original file line number Diff line number Diff line change
Expand Up @@ -2,10 +2,13 @@

# LQ Factorizations

struct LQ{T,S<:AbstractMatrix} <: Factorization{T}
struct LQ{T,S<:AbstractMatrix{T}} <: Factorization{T}
factors::S
τ::Vector{T}
LQ{T,S}(factors::AbstractMatrix{T}, τ::Vector{T}) where {T,S<:AbstractMatrix} = new(factors, τ)
function LQ{T,S}(factors, τ) where {T,S<:AbstractMatrix{T}}
@assert is_one_indexed(factors)
new(factors, τ)
end
end
LQ(factors::AbstractMatrix{T}, τ::Vector{T}) where {T} = LQ{T,typeof(factors)}(factors, τ)

Expand Down
7 changes: 5 additions & 2 deletions stdlib/LinearAlgebra/src/lu.jl
Original file line number Diff line number Diff line change
Expand Up @@ -3,11 +3,14 @@
####################
# LU Factorization #
####################
struct LU{T,S<:AbstractMatrix} <: Factorization{T}
struct LU{T,S<:AbstractMatrix{T}} <: Factorization{T}
factors::S
ipiv::Vector{BlasInt}
info::BlasInt
LU{T,S}(factors::AbstractMatrix{T}, ipiv::Vector{BlasInt}, info::BlasInt) where {T,S} = new(factors, ipiv, info)
function LU{T,S}(factors, ipiv, info) where {T,S<:AbstractMatrix{T}}
@assert is_one_indexed(factors)
new(factors, ipiv, info)
end
end
LU(factors::AbstractMatrix{T}, ipiv::Vector{BlasInt}, info::BlasInt) where {T} = LU{T,typeof(factors)}(factors, ipiv, info)

Expand Down
34 changes: 24 additions & 10 deletions stdlib/LinearAlgebra/src/qr.jl
Original file line number Diff line number Diff line change
Expand Up @@ -34,10 +34,13 @@ The object has two fields:
* `τ` is a vector of length `min(m,n)` containing the coefficients ``\tau_i``.
"""
struct QR{T,S<:AbstractMatrix} <: Factorization{T}
struct QR{T,S<:AbstractMatrix{T}} <: Factorization{T}
factors::S
τ::Vector{T}
QR{T,S}(factors::AbstractMatrix{T}, τ::Vector{T}) where {T,S<:AbstractMatrix} = new(factors, τ)
function QR{T,S}(factors, τ) where {T,S<:AbstractMatrix{T}}
@assert is_one_indexed(factors)
new(factors, τ)
end
end
QR(factors::AbstractMatrix{T}, τ::Vector{T}) where {T} = QR{T,typeof(factors)}(factors, τ)

Expand Down Expand Up @@ -94,10 +97,13 @@ The object has two fields:
[^Schreiber1989]: R Schreiber and C Van Loan, "A storage-efficient WY representation for products of Householder transformations", SIAM J Sci Stat Comput 10 (1989), 53-57. [doi:10.1137/0910005](https://doi.org/10.1137/0910005)
"""
struct QRCompactWY{S,M<:AbstractMatrix} <: Factorization{S}
struct QRCompactWY{S,M<:AbstractMatrix{S}} <: Factorization{S}
factors::M
T::Matrix{S}
QRCompactWY{S,M}(factors::AbstractMatrix{S}, T::AbstractMatrix{S}) where {S,M<:AbstractMatrix} = new(factors, T)
function QRCompactWY{S,M}(factors, T) where {S,M<:AbstractMatrix{S}}
@assert is_one_indexed(factors)
new(factors, T)
end
end
QRCompactWY(factors::AbstractMatrix{S}, T::AbstractMatrix{S}) where {S} = QRCompactWY{S,typeof(factors)}(factors, T)

Expand Down Expand Up @@ -139,12 +145,14 @@ The object has three fields:
* `jpvt` is an integer vector of length `n` corresponding to the permutation ``P``.
"""
struct QRPivoted{T,S<:AbstractMatrix} <: Factorization{T}
struct QRPivoted{T,S<:AbstractMatrix{T}} <: Factorization{T}
factors::S
τ::Vector{T}
jpvt::Vector{BlasInt}
QRPivoted{T,S}(factors::AbstractMatrix{T}, τ::Vector{T}, jpvt::Vector{BlasInt}) where {T,S<:AbstractMatrix} =
function QRPivoted{T,S}(factors, τ, jpvt) where {T,S<:AbstractMatrix{T}}
@assert is_one_indexed(factors)
new(factors, τ, jpvt)
end
end
QRPivoted(factors::AbstractMatrix{T}, τ::Vector{T}, jpvt::Vector{BlasInt}) where {T} =
QRPivoted{T,typeof(factors)}(factors, τ, jpvt)
Expand Down Expand Up @@ -435,10 +443,13 @@ abstract type AbstractQ{T} <: AbstractMatrix{T} end
The orthogonal/unitary ``Q`` matrix of a QR factorization stored in [`QR`](@ref) or
[`QRPivoted`](@ref) format.
"""
struct QRPackedQ{T,S<:AbstractMatrix} <: AbstractQ{T}
struct QRPackedQ{T,S<:AbstractMatrix{T}} <: AbstractQ{T}
factors::S
τ::Vector{T}
QRPackedQ{T,S}(factors::AbstractMatrix{T}, τ::Vector{T}) where {T,S<:AbstractMatrix} = new(factors, τ)
function QRPackedQ{T,S}(factors, τ) where {T,S<:AbstractMatrix{T}}
@assert is_one_indexed(factors)
new(factors, τ)
end
end
QRPackedQ(factors::AbstractMatrix{T}, τ::Vector{T}) where {T} = QRPackedQ{T,typeof(factors)}(factors, τ)

Expand All @@ -448,10 +459,13 @@ QRPackedQ(factors::AbstractMatrix{T}, τ::Vector{T}) where {T} = QRPackedQ{T,typ
The orthogonal/unitary ``Q`` matrix of a QR factorization stored in [`QRCompactWY`](@ref)
format.
"""
struct QRCompactWYQ{S, M<:AbstractMatrix} <: AbstractQ{S}
struct QRCompactWYQ{S, M<:AbstractMatrix{S}} <: AbstractQ{S}
factors::M
T::Matrix{S}
QRCompactWYQ{S,M}(factors::AbstractMatrix{S}, T::Matrix{S}) where {S,M<:AbstractMatrix} = new(factors, T)
function QRCompactWYQ{S,M}(factors, T) where {S,M<:AbstractMatrix{S}}
@assert is_one_indexed(factors)
new(factors, T)
end
end
QRCompactWYQ(factors::AbstractMatrix{S}, T::Matrix{S}) where {S} = QRCompactWYQ{S,typeof(factors)}(factors, T)

Expand Down
6 changes: 4 additions & 2 deletions stdlib/LinearAlgebra/src/svd.jl
Original file line number Diff line number Diff line change
@@ -1,12 +1,14 @@
# This file is a part of Julia. License is MIT: https://julialang.org/license

# Singular Value Decomposition
struct SVD{T,Tr,M<:AbstractArray} <: Factorization{T}
struct SVD{T,Tr,M<:AbstractArray{T}} <: Factorization{T}
U::M
S::Vector{Tr}
Vt::M
SVD{T,Tr,M}(U::AbstractArray{T}, S::Vector{Tr}, Vt::AbstractArray{T}) where {T,Tr,M} =
function SVD{T,Tr,M}(U, S, Vt) where {T,Tr,M<:AbstractArray{T}}
@assert is_one_indexed(U, S, Vt)
new(U, S, Vt)
end
end
SVD(U::AbstractArray{T}, S::Vector{Tr}, Vt::AbstractArray{T}) where {T,Tr} = SVD{T,Tr,typeof(U)}(U, S, Vt)

Expand Down
14 changes: 12 additions & 2 deletions stdlib/LinearAlgebra/src/symmetric.jl
Original file line number Diff line number Diff line change
@@ -1,9 +1,14 @@
# This file is a part of Julia. License is MIT: https://julialang.org/license

# Symmetric and Hermitian matrices
struct Symmetric{T,S<:AbstractMatrix{<:T}} <: AbstractMatrix{T}
struct Symmetric{T,S<:AbstractMatrix{T}} <: AbstractMatrix{T}
data::S
uplo::Char

function Symmetric{T,S}(data, uplo) where {T,S<:AbstractMatrix{T}}
@assert is_one_indexed(data)
new(data, uplo)
end
end
"""
Symmetric(A, uplo=:U)
Expand Down Expand Up @@ -71,9 +76,14 @@ function symmetric_type(::Type{T}) where {S, T<:AbstractMatrix{S}}
end
symmetric_type(::Type{T}) where {T<:Number} = T

struct Hermitian{T,S<:AbstractMatrix{<:T}} <: AbstractMatrix{T}
struct Hermitian{T,S<:AbstractMatrix{T}} <: AbstractMatrix{T}
data::S
uplo::Char

function Hermitian{T,S}(data, uplo) where {T,S<:AbstractMatrix{T}}
@assert is_one_indexed(data)
new(data, uplo)
end
end
"""
Hermitian(A, uplo=:U)
Expand Down
9 changes: 7 additions & 2 deletions stdlib/LinearAlgebra/src/triangular.jl
Original file line number Diff line number Diff line change
Expand Up @@ -9,13 +9,18 @@ abstract type AbstractTriangular{T,S<:AbstractMatrix} <: AbstractMatrix{T} end
for t in (:LowerTriangular, :UnitLowerTriangular, :UpperTriangular,
:UnitUpperTriangular)
@eval begin
struct $t{T,S<:AbstractMatrix} <: AbstractTriangular{T,S}
struct $t{T,S<:AbstractMatrix{T}} <: AbstractTriangular{T,S}
data::S

function $t{T,S}(data) where {T,S<:AbstractMatrix{T}}
@assert is_one_indexed(data)
checksquare(data)
new(data)
end
end
$t(A::$t) = A
$t{T}(A::$t{T}) where {T} = A
function $t(A::AbstractMatrix)
checksquare(A)
return $t{eltype(A), typeof(A)}(A)
end

Expand Down
6 changes: 4 additions & 2 deletions stdlib/LinearAlgebra/src/tridiag.jl
Original file line number Diff line number Diff line change
Expand Up @@ -6,11 +6,12 @@
struct SymTridiagonal{T,V<:AbstractVector{T}} <: AbstractMatrix{T}
dv::V # diagonal
ev::V # subdiagonal
function SymTridiagonal{T}(dv::V, ev::V) where {T,V<:AbstractVector{T}}
function SymTridiagonal{T,V}(dv, ev) where {T,V<:AbstractVector{T}}
@assert is_one_indexed(dv, ev)
if !(length(dv) - 1 <= length(ev) <= length(dv))
throw(DimensionMismatch("subdiagonal has wrong length. Has length $(length(ev)), but should be either $(length(dv) - 1) or $(length(dv))."))
end
new{T,V}(dv,ev)
new(dv,ev)
end
end

Expand Down Expand Up @@ -46,6 +47,7 @@ julia> SymTridiagonal(dv, ev)
```
"""
SymTridiagonal(dv::V, ev::V) where {T,V<:AbstractVector{T}} = SymTridiagonal{T}(dv, ev)
SymTridiagonal{T}(dv::V, ev::V) where {T,V<:AbstractVector{T}} = SymTridiagonal{T,V}(dv, ev)

"""
SymTridiagonal(A::AbstractMatrix)
Expand Down

0 comments on commit da19abe

Please sign in to comment.