Skip to content

Commit 72ca30a

Browse files
author
Andy Ferris
committed
Added ConjArray wrapper type for conjugate views
By default, this is used only in conjugation with `RowVector`, so that both `transpose(vec)` and `ctranspose(vec)` both return views.
1 parent 4186749 commit 72ca30a

File tree

8 files changed

+80
-4
lines changed

8 files changed

+80
-4
lines changed

NEWS.md

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -185,6 +185,10 @@ Library improvements
185185

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

188+
* Introduced a wrapper type for lazy complex conjugation of arrays, `ConjArray`.
189+
Currently, it is used by default for the new `RowVector` type only, and
190+
enforces that both `transpose(vec)` and `ctranspose(vec)` are views not copies.
191+
188192
Compiler/Runtime improvements
189193
-----------------------------
190194

base/exports.jl

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -49,6 +49,7 @@ export
4949
Complex128,
5050
Complex64,
5151
Complex32,
52+
ConjArray,
5253
DenseMatrix,
5354
DenseVecOrMat,
5455
DenseVector,

base/linalg/conjarray.jl

Lines changed: 45 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,45 @@
1+
# This file is a part of Julia. License is MIT: http://julialang.org/license
2+
3+
"""
4+
ConjArray(array)
5+
6+
A lazy-view wrapper of an `AbstractArray`, taking the elementwise complex
7+
conjugate. This type is usually constructed (and unwrapped) via the `conj()`
8+
function (or related `ctranspose()`), but currently this is the default behavior
9+
for `RowVector` only.
10+
"""
11+
immutable ConjArray{T, N, A <: AbstractArray} <: AbstractArray{T, N}
12+
parent::A
13+
end
14+
15+
@inline ConjArray{T,N}(a::AbstractArray{T,N}) = ConjArray{conj_type(T), N, typeof(a)}(a)
16+
17+
# This type can cause the element type to change under conjugation - e.g. an array of complex arrays.
18+
@inline conj_type(x) = conj_type(typeof(x))
19+
@inline conj_type{T}(::Type{T}) = promote_op(conj, T)
20+
21+
@inline parent(c::ConjArray) = c.parent
22+
@inline parent_type(c::ConjArray) = parent_type(typeof(c))
23+
@inline parent_type{T,N,A}(::Type{ConjArray{T,N,A}}) = A
24+
25+
@inline size(a::ConjArray) = size(a.parent)
26+
linearindexing{CA <: ConjArray}(::Union{CA,Type{CA}}) = linearindexing(parent_type(CA))
27+
28+
@propagate_inbounds getindex{T,N}(a::ConjArray{T,N}, i::Int) = conj(getindex(a.parent, i))
29+
@propagate_inbounds getindex{T,N}(a::ConjArray{T,N}, i::Vararg{Int,N}) = conj(getindex(a.parent, i...))
30+
@propagate_inbounds setindex!{T,N}(a::ConjArray{T,N}, v, i::Int) = setindex!(a.parent, conj(v), i)
31+
@propagate_inbounds setindex!{T,N}(a::ConjArray{T,N}, v, i::Vararg{Int,N}) = setindex!(a.parent, conj(v), i...)
32+
33+
# Currently, this is default behavior for RowVector only
34+
"""
35+
conj(rowvector)
36+
37+
Returns a `ConjArray` lazy view of the input, where each element is conjugated.
38+
"""
39+
@inline conj(rv::RowVector) = RowVector(_conj(parent(rv)))
40+
@inline conj(a::ConjArray) = parent(a)
41+
42+
@inline _conj(a::AbstractArray) = ConjArray(a)
43+
@inline _conj{T<:Real}(a::AbstractArray{T}) = a
44+
@inline _conj(a::ConjArray) = parent(a)
45+
@inline _conj{T<:Real}(a::ConjArray{T}) = parent(a)

base/linalg/linalg.jl

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -22,6 +22,7 @@ export
2222

2323
# Types
2424
RowVector,
25+
ConjArray,
2526
SymTridiagonal,
2627
Tridiagonal,
2728
Bidiagonal,
@@ -239,6 +240,7 @@ copy_oftype{T,N,S}(A::AbstractArray{T,N}, ::Type{S}) = convert(AbstractArray{S,N
239240

240241
include("transpose.jl")
241242
include("rowvector.jl")
243+
include("conjarray.jl")
242244

243245
include("exceptions.jl")
244246
include("generic.jl")

base/linalg/rowvector.jl

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -67,18 +67,18 @@ julia> transpose(v)
6767
```
6868
"""
6969
@inline transpose(vec::AbstractVector) = RowVector(vec)
70-
@inline ctranspose{T}(vec::AbstractVector{T}) = RowVector(conj(vec))
70+
@inline ctranspose{T}(vec::AbstractVector{T}) = RowVector(_conj(vec))
7171
@inline ctranspose{T<:Real}(vec::AbstractVector{T}) = RowVector(vec)
7272

7373
@inline transpose(rowvec::RowVector) = rowvec.vec
74-
@inline ctranspose{T}(rowvec::RowVector{T}) = conj(rowvec.vec)
74+
@inline ctranspose{T}(rowvec::RowVector{T}) = _conj(rowvec.vec)
7575
@inline ctranspose{T<:Real}(rowvec::RowVector{T}) = rowvec.vec
7676

7777
parent(rowvec::RowVector) = rowvec.vec
7878

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

8484
# AbstractArray interface

doc/src/stdlib/linalg.md

Lines changed: 1 addition & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -100,6 +100,7 @@ Base.LinAlg.istriu
100100
Base.LinAlg.isdiag
101101
Base.LinAlg.ishermitian
102102
Base.LinAlg.RowVector
103+
Base.LinAlg.ConjArray
103104
Base.transpose
104105
Base.transpose!
105106
Base.ctranspose

test/choosetests.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -130,7 +130,7 @@ function choosetests(choices = [])
130130
"linalg/diagonal", "linalg/pinv", "linalg/givens",
131131
"linalg/cholesky", "linalg/lu", "linalg/symmetric",
132132
"linalg/generic", "linalg/uniformscaling", "linalg/lq",
133-
"linalg/hessenberg", "linalg/rowvector"]
133+
"linalg/hessenberg", "linalg/rowvector", "linalg/conjarray"]
134134
if Base.USE_GPL_LIBS
135135
push!(linalgtests, "linalg/arnoldi")
136136
end

test/linalg/conjarray.jl

Lines changed: 23 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,23 @@
1+
# This file is a part of Julia. License is MIT: http://julialang.org/license
2+
3+
@testset "RowVector" begin
4+
5+
@testset "Core" begin
6+
m = [1+im 2; 2 4-im]
7+
cm = ConjArray(m)
8+
@test cm[1,1] == 1-im
9+
@test trace(cm*m) == 27
10+
11+
v = [[1+im], [1-im]]
12+
cv = ConjArray(v)
13+
@test cv[1] == [1-im]
14+
end
15+
16+
@testset "RowVector conjugates" begin
17+
v = [1+im, 1-im]
18+
rv = v'
19+
@test (parent(rv) isa ConjArray)
20+
@test rv' === v
21+
end
22+
23+
end

0 commit comments

Comments
 (0)