Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Extend FTPlans to Arrays and higher dimensions #252

Merged
Show file tree
Hide file tree
Changes from 6 commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions src/FastTransforms.jl
Original file line number Diff line number Diff line change
Expand Up @@ -128,6 +128,7 @@ for f in (:jac2jac,
@eval $f(x::AbstractArray, y...; z...) = $lib_f(x, y...; z...)
end

include("arrays.jl")
# following use Toeplitz-Hankel to avoid expensive plans
# for f in (:leg2cheb, :cheb2leg, :ultra2ultra)
# th_f = Symbol("th_", f)
Expand Down
90 changes: 90 additions & 0 deletions src/arrays.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,90 @@
struct ArrayPlan{T, FF<:FTPlan{<:T}, Szs<:Tuple, Dims<:Tuple{<:Int}} <: Plan{T}
F::FF
szs::Szs
dims::Dims
end
size(P::ArrayPlan) = P.szs
size(P::ArrayPlan, k::Int) = P.szs[k]
size(P::ArrayPlan, k...) = P.szs[[k...]]

function ArrayPlan(F::FTPlan{<:T}, c::AbstractArray{T}, dims::Tuple{<:Int}=(1,)) where T
szs = size(c)
@assert F.n == szs[dims[1]]
ArrayPlan(F, size(c), dims)
end

function *(P::ArrayPlan, f::AbstractArray)
F, dims, szs = P.F, P.dims, P.szs
@assert length(dims) == 1
@assert szs == size(f)
d = first(dims)

perm = (d, ntuple(i-> i + (i >= d), ndims(f) -1)...)
fp = permutedims(f, perm)

fr = reshape(fp, size(fp,1), :)

permutedims(reshape(F*fr, size(fp)...), invperm(perm))
end

function \(P::ArrayPlan, f::AbstractArray)
F, dims, szs = P.F, P.dims, P.szs
@assert length(dims) == 1
@assert szs == size(f)
d = first(dims)

perm = (d, ntuple(i-> i + (i >= d), ndims(f) -1)...)
fp = permutedims(f, perm)

fr = reshape(fp, size(fp,1), :)

permutedims(reshape(F\fr, size(fp)...), invperm(perm))
end

struct NDimsPlan{T, FF<:ArrayPlan{<:T}, Szs<:Tuple, Dims<:Tuple} <: Plan{T}
F::FF
szs::Szs
dims::Dims
function NDimsPlan(F, szs, dims)
if length(Set(szs[[dims...]])) > 1
error("Different size in dims axes not yet implemented in N-dimensional transform.")
end
new{eltype(F), typeof(F), typeof(szs), typeof(dims)}(F, szs, dims)
end
end

size(P::NDimsPlan) = P.szs
size(P::NDimsPlan, k::Int) = P.szs[k]
size(P::NDimsPlan, k...) = P.szs[[k...]]
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Why the array construction in the indexing? Shouldn't this just be P.szs[k...]?

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

I've probably misunderstood what k is. I copied the size(P, k...) structure from size(P::JacobiTransformPlan, k...) = size(P.chebtransform, k...) and assumed it is a tuple. In which case P.szs[k...] does not work.

e.g.

julia> szs = (10,20,30);
julia> szs[[(1,2)...]]
(10, 20)
julia> szs[(1,2)...]
ERROR: MethodError: no method matching getindex(::Tuple{Int64, Int64, Int64}, ::Int64, ::Int64)

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Perhaps I've misunderstood this. Are you trying to obtain multiple indices in one go?

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

That's the intention! Although I don't need that functionality in this PR.

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Perhaps best to remove it then, as it's not the conventional usage of size. One may broadcast size over indices to obtain the sizes along multiple dimensions at once.

Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Note sizes of plans and arrays behave differently

compare w FFT

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

judging by this line here https://github.com/JuliaMath/FFTW.jl/blob/f888022d7a1ff78491abf8f33f1055cc52a68f0a/src/fft.jl#L254 I think I should only really keep size(P::NDimsPlan) = P.szs

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

...or maybe I should keep size(P::NDimsPlan, k::Int) = P.szs[k]?


function NDimsPlan(F::FTPlan, szs::Tuple, dims::Tuple)
NDimsPlan(ArrayPlan(F, szs, (first(dims),)), szs, dims)
end

function *(P::NDimsPlan, f::AbstractArray)
F, dims = P.F, P.dims
@assert size(P) == size(f)
g = copy(f)
t = 1:ndims(g)
d1 = dims[1]
for d in dims
perm = ntuple(k -> k == d1 ? t[d] : k == d ? t[d1] : t[k], ndims(g))
gp = permutedims(g, perm)
g = permutedims(F*gp, invperm(perm))
end
return g
end

function \(P::NDimsPlan, f::AbstractArray)
F, dims = P.F, P.dims
@assert size(P) == size(f)
g = copy(f)
t = 1:ndims(g)
d1 = dims[1]
for d in dims
perm = ntuple(k -> k == d1 ? t[d] : k == d ? t[d1] : t[k], ndims(g))
gp = permutedims(g, perm)
g = permutedims(F\gp, invperm(perm))
end
return g
end
68 changes: 68 additions & 0 deletions test/arraystests.jl
Original file line number Diff line number Diff line change
@@ -0,0 +1,68 @@
using FastTransforms, Test
import FastTransforms: ArrayPlan, NDimsPlan

@testset "Array transform" begin
@testset "ArrayPlan" begin
c = randn(5,20,10)
F = plan_cheb2leg(c)
FT = ArrayPlan(F, c)

@test size(FT) == size(c)
@test size(FT,1) == size(c,1)
@test size(FT,1,2) == (size(c,1), size(c,2))

f = similar(c);
for k in axes(c,3)
f[:,:,k] = (F*c[:,:,k])
end
@test f ≈ FT*c
@test c ≈ FT\f

F = plan_cheb2leg(Vector{Float64}(axes(c,2)))
FT = ArrayPlan(F, c, (2,))
for k in axes(c,3)
f[:,:,k] = (F*c[:,:,k]')'
end
@test f ≈ FT*c
@test c ≈ FT\f
end

@testset "NDimsPlan" begin
c = randn(20,10,20)
@test_throws ErrorException("Different size in dims axes not yet implemented in N-dimensional transform.") NDimsPlan(ArrayPlan(plan_cheb2leg(c), c), size(c), (1,2))

c = randn(5,20)
F = plan_cheb2leg(c)
FT = ArrayPlan(F, c)
P = NDimsPlan(F, size(c), (1,))
@test F*c ≈ FT*c ≈ P*c

c = randn(20,20,5);
F = plan_cheb2leg(c)
FT = ArrayPlan(F, c)
P = NDimsPlan(FT, size(c), (1,2))

@test size(P) == size(c)
@test size(P,1) == size(c,1)
@test size(P,1,2) == (size(c,1), size(c,2))

f = similar(c);
for k in axes(f,3)
f[:,:,k] = (F*(F*c[:,:,k])')'
end
@test f ≈ P*c
@test c ≈ P\f

c = randn(5,10,10,60)
F = plan_cheb2leg(randn(10))
P = NDimsPlan(F, size(c), (2,3))
f = similar(c)
for i in axes(f,1), j in axes(f,4)
f[i,:,:,j] = (F*(F*c[i,:,:,j])')'
end
@test f ≈ P*c
@test c ≈ P\f
end
end


1 change: 1 addition & 0 deletions test/runtests.jl
Original file line number Diff line number Diff line change
Expand Up @@ -12,3 +12,4 @@ include("clenshawtests.jl")
include("toeplitzplanstests.jl")
include("toeplitzhankeltests.jl")
include("symmetrictoeplitzplushankeltests.jl")
include("arraystests.jl")
Loading