Skip to content

Commit a40b3dc

Browse files
authored
Change clenshaw!/horner! to mutate first argument (#255)
* Change clenshaw!/horner! to mutate first argument * Update libfasttransformstests.jl * Update Project.toml * Update make.jl * Update sphere.jl * Update make.jl * Move semiclassical.jl to ApproxFunExamples
1 parent 5260921 commit a40b3dc

8 files changed

+20
-140
lines changed

Project.toml

+3-4
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,7 @@
11
name = "FastTransforms"
22
uuid = "057dd010-8810-581a-b7be-e3fc3b93f78c"
3-
version = "0.16.8"
3+
version = "0.17"
4+
45

56
[deps]
67
AbstractFFTs = "621f4979-c628-5d54-868e-fcf4e3e8185c"
@@ -15,7 +16,6 @@ LazyArrays = "5078a376-72f3-5289-bfd5-ec5146d43c02"
1516
Libdl = "8f399da3-3557-5675-b5ff-fb832c97cbdb"
1617
LinearAlgebra = "37e2e46d-f89d-539d-b4ee-838fcccc9c8e"
1718
RecurrenceRelationships = "807425ed-42ea-44d6-a357-6771516d7b2c"
18-
Reexport = "189a3867-3050-52da-a836-e630ba90ab69"
1919
SpecialFunctions = "276daf66-3868-5448-9aa4-cd146d93841b"
2020
ToeplitzMatrices = "c751599d-da0a-543b-9d20-d0a503d91d24"
2121

@@ -29,8 +29,7 @@ FastTransforms_jll = "0.6.2"
2929
FillArrays = "0.9, 0.10, 0.11, 0.12, 0.13, 1"
3030
GenericFFT = "0.1"
3131
LazyArrays = "2.2"
32-
RecurrenceRelationships = "0.1"
33-
Reexport = "0.2, 1.0"
32+
RecurrenceRelationships = "0.2"
3433
SpecialFunctions = "0.10, 1, 2"
3534
ToeplitzMatrices = "0.7.1, 0.8"
3635
julia = "1.7"

docs/Project.toml

-1
Original file line numberDiff line numberDiff line change
@@ -1,5 +1,4 @@
11
[deps]
2-
ApproxFun = "28f2ccd6-bb30-5033-b560-165f7b14dc2f"
32
Documenter = "e30172f5-a6a5-5a46-863b-614d45cd2de4"
43
FastTransforms = "057dd010-8810-581a-b7be-e3fc3b93f78c"
54
LaTeXStrings = "b964fa9f-0449-5b57-a5c2-d3ea65f4040f"

docs/make.jl

-2
Original file line numberDiff line numberDiff line change
@@ -13,7 +13,6 @@ examples = [
1313
"halfrange.jl",
1414
"nonlocaldiffusion.jl",
1515
"padua.jl",
16-
"semiclassical.jl",
1716
"sphere.jl",
1817
"spinweighted.jl",
1918
"subspaceangles.jl",
@@ -48,7 +47,6 @@ makedocs(
4847
"generated/halfrange.md",
4948
"generated/nonlocaldiffusion.md",
5049
"generated/padua.md",
51-
"generated/semiclassical.md",
5250
"generated/sphere.md",
5351
"generated/spinweighted.md",
5452
"generated/subspaceangles.md",

examples/semiclassical.jl

-113
This file was deleted.

examples/sphere.jl

+3-3
Original file line numberDiff line numberDiff line change
@@ -61,7 +61,7 @@ C = [k/(k+1) for k in 0:N]
6161
c = zeros(N); c[N] = 1
6262
pts = vec([z(θ, φ)y for θ in θ, φ in φ])
6363
phi0 = ones(N*M)
64-
F = reshape(FastTransforms.clenshaw!(c, A, B, C, pts, phi0, zeros(N*M)), N, M)
64+
F = reshape(FastTransforms.clenshaw!(zeros(N*M), c, A, B, C, pts, phi0), N, M)
6565

6666
# We superpose a surface plot of $f$ on top of the grid:
6767
X = [sinpi(θ)*cospi(φ) for θ in θ, φ in φ]
@@ -91,7 +91,7 @@ U = threshold!(P\V, 400*eps())
9191
nrm1 = norm(U)
9292

9393
# Similarly, on the tensor product grid, our function samples are:
94-
Pnxy = FastTransforms.clenshaw!(c, A, B, C, [xy], [1.0], [0.0])[1]
94+
Pnxy = FastTransforms.clenshaw!([0.0], c, A, B, C, [xy], [1.0])[1]
9595
F = [(F[n, m] - Pnxy)/(z(θ[n], φ[m])y - xy) for n in 1:N, m in 1:M]
9696

9797
# We superpose a surface plot of $f$ on top of the grid:
@@ -108,7 +108,7 @@ U = threshold!(P\V, 400*eps())
108108

109109
# Finally, the Legendre polynomial $P_n(z\cdot x)$ is aligned with the grid:
110110
pts = vec([z(θ, φ)x for θ in θ, φ in φ])
111-
F = reshape(FastTransforms.clenshaw!(c, A, B, C, pts, phi0, zeros(N*M)), N, M)
111+
F = reshape(FastTransforms.clenshaw!(zeros(N*M), c, A, B, C, pts, phi0), N, M)
112112

113113
# We superpose a surface plot of $f$ on top of the grid:
114114
scatter3d(vec(X), vec(Y), vec(Z); markersize=1.25, markercolor=:violetred)

src/FastTransforms.jl

+5-8
Original file line numberDiff line numberDiff line change
@@ -1,11 +1,11 @@
11
module FastTransforms
22

33
using ArrayLayouts, BandedMatrices, FastGaussQuadrature, FillArrays, LazyArrays, LinearAlgebra,
4-
Reexport, SpecialFunctions, ToeplitzMatrices, RecurrenceRelationships
4+
SpecialFunctions, ToeplitzMatrices, RecurrenceRelationships
55

6-
@reexport using AbstractFFTs
7-
@reexport using FFTW
8-
@reexport using GenericFFT
6+
using AbstractFFTs
7+
using FFTW
8+
using GenericFFT
99

1010
import Base: convert, unsafe_convert, eltype, ndims, adjoint, transpose, show,
1111
*, \, inv, length, size, view, getindex, tail, OneTo
@@ -34,11 +34,8 @@ import LinearAlgebra: cholesky, issymmetric, isposdef, mul!, lmul!, ldiv!
3434

3535
import GenericFFT: interlace # imported in downstream packages
3636

37-
import RecurrenceRelationships: clenshaw!, check_clenshaw_recurrences
37+
import RecurrenceRelationships: check_clenshaw_recurrences
3838

39-
const _forwardrecurrence! = RecurrenceRelationships.forwardrecurrence!
40-
const _clenshaw_next = RecurrenceRelationships.clenshaw_next
41-
const _forwardrecurrence_next = RecurrenceRelationships.forwardrecurrence_next
4239

4340
export leg2cheb, cheb2leg, ultra2ultra, jac2jac,
4441
lag2lag, jac2ultra, ultra2jac, jac2cheb,

src/libfasttransforms.jl

+6-6
Original file line numberDiff line numberDiff line change
@@ -49,13 +49,13 @@ function renew!(x::AbstractArray{BigFloat})
4949
return x
5050
end
5151

52-
function horner!(c::StridedVector{Float64}, x::Vector{Float64}, f::Vector{Float64})
52+
function horner!(f::Vector{Float64}, c::StridedVector{Float64}, x::Vector{Float64})
5353
@assert length(x) == length(f)
5454
ccall((:ft_horner, libfasttransforms), Cvoid, (Cint, Ptr{Float64}, Cint, Cint, Ptr{Float64}, Ptr{Float64}), length(c), c, stride(c, 1), length(x), x, f)
5555
f
5656
end
5757

58-
function horner!(c::StridedVector{Float32}, x::Vector{Float32}, f::Vector{Float32})
58+
function horner!(f::Vector{Float32}, c::StridedVector{Float32}, x::Vector{Float32})
5959
@assert length(x) == length(f)
6060
ccall((:ft_hornerf, libfasttransforms), Cvoid, (Cint, Ptr{Float32}, Cint, Cint, Ptr{Float32}, Ptr{Float32}), length(c), c, stride(c, 1), length(x), x, f)
6161
f
@@ -69,27 +69,27 @@ function check_clenshaw_points(x, f)
6969
length(x) == length(f) || throw(ArgumentError("Dimensions must match"))
7070
end
7171

72-
function clenshaw!(c::StridedVector{Float64}, x::Vector{Float64}, f::Vector{Float64})
72+
function clenshaw!(f::Vector{Float64}, c::StridedVector{Float64}, x::Vector{Float64})
7373
@boundscheck check_clenshaw_points(x, f)
7474
ccall((:ft_clenshaw, libfasttransforms), Cvoid, (Cint, Ptr{Float64}, Cint, Cint, Ptr{Float64}, Ptr{Float64}), length(c), c, stride(c, 1), length(x), x, f)
7575
f
7676
end
7777

78-
function clenshaw!(c::StridedVector{Float32}, x::Vector{Float32}, f::Vector{Float32})
78+
function clenshaw!(f::Vector{Float32}, c::StridedVector{Float32}, x::Vector{Float32})
7979
@boundscheck check_clenshaw_points(x, f)
8080
ccall((:ft_clenshawf, libfasttransforms), Cvoid, (Cint, Ptr{Float32}, Cint, Cint, Ptr{Float32}, Ptr{Float32}), length(c), c, stride(c, 1), length(x), x, f)
8181
f
8282
end
8383

84-
function clenshaw!(c::StridedVector{Float64}, A::Vector{Float64}, B::Vector{Float64}, C::Vector{Float64}, x::Vector{Float64}, ϕ₀::Vector{Float64}, f::Vector{Float64})
84+
function clenshaw!(f::Vector{Float64}, c::StridedVector{Float64}, A::Vector{Float64}, B::Vector{Float64}, C::Vector{Float64}, x::Vector{Float64}, ϕ₀::Vector{Float64})
8585
N = length(c)
8686
@boundscheck check_clenshaw_recurrences(N, A, B, C)
8787
@boundscheck check_clenshaw_points(x, ϕ₀, f)
8888
ccall((:ft_orthogonal_polynomial_clenshaw, libfasttransforms), Cvoid, (Cint, Ptr{Float64}, Cint, Ptr{Float64}, Ptr{Float64}, Ptr{Float64}, Cint, Ptr{Float64}, Ptr{Float64}, Ptr{Float64}), N, c, stride(c, 1), A, B, C, length(x), x, ϕ₀, f)
8989
f
9090
end
9191

92-
function clenshaw!(c::StridedVector{Float32}, A::Vector{Float32}, B::Vector{Float32}, C::Vector{Float32}, x::Vector{Float32}, ϕ₀::Vector{Float32}, f::Vector{Float32})
92+
function clenshaw!(f::Vector{Float32}, c::StridedVector{Float32}, A::Vector{Float32}, B::Vector{Float32}, C::Vector{Float32}, x::Vector{Float32}, ϕ₀::Vector{Float32})
9393
N = length(c)
9494
@boundscheck check_clenshaw_recurrences(N, A, B, C)
9595
@boundscheck check_clenshaw_points(x, ϕ₀, f)

test/libfasttransformstests.jl

+3-3
Original file line numberDiff line numberDiff line change
@@ -8,18 +8,18 @@ FastTransforms.ft_set_num_threads(ceil(Int, Base.Sys.CPU_THREADS/2))
88
c = one(T) ./ (1:n)
99
x = collect(-1 .+ 2*(0:n-1)/T(n))
1010
f = similar(x)
11-
@test FastTransforms.horner!(c, x, f) == f
11+
@test FastTransforms.horner!(f, c, x) == f
1212
fd = T[sum(c[k]*x^(k-1) for k in 1:length(c)) for x in x]
1313
@test f fd
14-
@test FastTransforms.clenshaw!(c, x, f) == f
14+
@test FastTransforms.clenshaw!(f, c, x) == f
1515
fd = T[sum(c[k]*cos((k-1)*acos(x)) for k in 1:length(c)) for x in x]
1616
@test f fd
1717
A = T[(2k+one(T))/(k+one(T)) for k in 0:length(c)-1]
1818
B = T[zero(T) for k in 0:length(c)-1]
1919
C = T[k/(k+one(T)) for k in 0:length(c)]
2020
phi0 = ones(T, length(x))
2121
c = FastTransforms.lib_cheb2leg(c)
22-
@test FastTransforms.clenshaw!(c, A, B, C, x, phi0, f) == f
22+
@test FastTransforms.clenshaw!(f, c, A, B, C, x, phi0) == f
2323
@test f fd
2424
end
2525

0 commit comments

Comments
 (0)