Skip to content

Commit

Permalink
Added ConjArray wrapper type for conjugate views
Browse files Browse the repository at this point in the history
By default, this is used only in conjugation with `RowVector`, so that
both `transpose(vec)` and `ctranspose(vec)` both return views.
  • Loading branch information
Andy Ferris committed Jan 15, 2017
1 parent 4186749 commit 72ca30a
Show file tree
Hide file tree
Showing 8 changed files with 80 additions and 4 deletions.
4 changes: 4 additions & 0 deletions NEWS.md
Original file line number Diff line number Diff line change
Expand Up @@ -185,6 +185,10 @@ Library improvements

* `notify` now returns a count of tasks woken up ([#19841]).

* Introduced a wrapper type for lazy complex conjugation of arrays, `ConjArray`.
Currently, it is used by default for the new `RowVector` type only, and
enforces that both `transpose(vec)` and `ctranspose(vec)` are views not copies.

Compiler/Runtime improvements
-----------------------------

Expand Down
1 change: 1 addition & 0 deletions base/exports.jl
Original file line number Diff line number Diff line change
Expand Up @@ -49,6 +49,7 @@ export
Complex128,
Complex64,
Complex32,
ConjArray,
DenseMatrix,
DenseVecOrMat,
DenseVector,
Expand Down
45 changes: 45 additions & 0 deletions base/linalg/conjarray.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,45 @@
# This file is a part of Julia. License is MIT: http://julialang.org/license

"""
ConjArray(array)
A lazy-view wrapper of an `AbstractArray`, taking the elementwise complex
conjugate. This type is usually constructed (and unwrapped) via the `conj()`
function (or related `ctranspose()`), but currently this is the default behavior
for `RowVector` only.
"""
immutable ConjArray{T, N, A <: AbstractArray} <: AbstractArray{T, N}
parent::A
end

@inline ConjArray{T,N}(a::AbstractArray{T,N}) = ConjArray{conj_type(T), N, typeof(a)}(a)

# This type can cause the element type to change under conjugation - e.g. an array of complex arrays.
@inline conj_type(x) = conj_type(typeof(x))
@inline conj_type{T}(::Type{T}) = promote_op(conj, T)

@inline parent(c::ConjArray) = c.parent
@inline parent_type(c::ConjArray) = parent_type(typeof(c))
@inline parent_type{T,N,A}(::Type{ConjArray{T,N,A}}) = A

@inline size(a::ConjArray) = size(a.parent)
linearindexing{CA <: ConjArray}(::Union{CA,Type{CA}}) = linearindexing(parent_type(CA))

@propagate_inbounds getindex{T,N}(a::ConjArray{T,N}, i::Int) = conj(getindex(a.parent, i))
@propagate_inbounds getindex{T,N}(a::ConjArray{T,N}, i::Vararg{Int,N}) = conj(getindex(a.parent, i...))
@propagate_inbounds setindex!{T,N}(a::ConjArray{T,N}, v, i::Int) = setindex!(a.parent, conj(v), i)
@propagate_inbounds setindex!{T,N}(a::ConjArray{T,N}, v, i::Vararg{Int,N}) = setindex!(a.parent, conj(v), i...)

# Currently, this is default behavior for RowVector only
"""
conj(rowvector)
Returns a `ConjArray` lazy view of the input, where each element is conjugated.
"""
@inline conj(rv::RowVector) = RowVector(_conj(parent(rv)))
@inline conj(a::ConjArray) = parent(a)

@inline _conj(a::AbstractArray) = ConjArray(a)
@inline _conj{T<:Real}(a::AbstractArray{T}) = a
@inline _conj(a::ConjArray) = parent(a)
@inline _conj{T<:Real}(a::ConjArray{T}) = parent(a)
2 changes: 2 additions & 0 deletions base/linalg/linalg.jl
Original file line number Diff line number Diff line change
Expand Up @@ -22,6 +22,7 @@ export

# Types
RowVector,
ConjArray,
SymTridiagonal,
Tridiagonal,
Bidiagonal,
Expand Down Expand Up @@ -239,6 +240,7 @@ copy_oftype{T,N,S}(A::AbstractArray{T,N}, ::Type{S}) = convert(AbstractArray{S,N

include("transpose.jl")
include("rowvector.jl")
include("conjarray.jl")

include("exceptions.jl")
include("generic.jl")
Expand Down
6 changes: 3 additions & 3 deletions base/linalg/rowvector.jl
Original file line number Diff line number Diff line change
Expand Up @@ -67,18 +67,18 @@ julia> transpose(v)
```
"""
@inline transpose(vec::AbstractVector) = RowVector(vec)
@inline ctranspose{T}(vec::AbstractVector{T}) = RowVector(conj(vec))
@inline ctranspose{T}(vec::AbstractVector{T}) = RowVector(_conj(vec))
@inline ctranspose{T<:Real}(vec::AbstractVector{T}) = RowVector(vec)

@inline transpose(rowvec::RowVector) = rowvec.vec
@inline ctranspose{T}(rowvec::RowVector{T}) = conj(rowvec.vec)
@inline ctranspose{T}(rowvec::RowVector{T}) = _conj(rowvec.vec)
@inline ctranspose{T<:Real}(rowvec::RowVector{T}) = rowvec.vec

parent(rowvec::RowVector) = rowvec.vec

# Strictly, these are unnecessary but will make things stabler if we introduce
# a "view" for conj(::AbstractArray)
@inline conj(rowvec::RowVector) = RowVector(conj(rowvec.vec))
@inline conj(rowvec::RowVector) = RowVector(_conj(rowvec.vec))
@inline conj{T<:Real}(rowvec::RowVector{T}) = rowvec

# AbstractArray interface
Expand Down
1 change: 1 addition & 0 deletions doc/src/stdlib/linalg.md
Original file line number Diff line number Diff line change
Expand Up @@ -100,6 +100,7 @@ Base.LinAlg.istriu
Base.LinAlg.isdiag
Base.LinAlg.ishermitian
Base.LinAlg.RowVector
Base.LinAlg.ConjArray
Base.transpose
Base.transpose!
Base.ctranspose
Expand Down
2 changes: 1 addition & 1 deletion test/choosetests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -130,7 +130,7 @@ function choosetests(choices = [])
"linalg/diagonal", "linalg/pinv", "linalg/givens",
"linalg/cholesky", "linalg/lu", "linalg/symmetric",
"linalg/generic", "linalg/uniformscaling", "linalg/lq",
"linalg/hessenberg", "linalg/rowvector"]
"linalg/hessenberg", "linalg/rowvector", "linalg/conjarray"]
if Base.USE_GPL_LIBS
push!(linalgtests, "linalg/arnoldi")
end
Expand Down
23 changes: 23 additions & 0 deletions test/linalg/conjarray.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,23 @@
# This file is a part of Julia. License is MIT: http://julialang.org/license

@testset "RowVector" begin

@testset "Core" begin
m = [1+im 2; 2 4-im]
cm = ConjArray(m)
@test cm[1,1] == 1-im
@test trace(cm*m) == 27

v = [[1+im], [1-im]]
cv = ConjArray(v)
@test cv[1] == [1-im]
end

@testset "RowVector conjugates" begin
v = [1+im, 1-im]
rv = v'
@test (parent(rv) isa ConjArray)
@test rv' === v
end

end

0 comments on commit 72ca30a

Please sign in to comment.