Skip to content

Commit 7ca43ce

Browse files
add Symmetric(Banded)ToeplitzPlusHankel
1 parent 9f135a8 commit 7ca43ce

5 files changed

+185
-4
lines changed

Project.toml

+3-1
Original file line numberDiff line numberDiff line change
@@ -1,9 +1,10 @@
11
name = "FastTransforms"
22
uuid = "057dd010-8810-581a-b7be-e3fc3b93f78c"
3-
version = "0.15.16"
3+
version = "0.16.0"
44

55
[deps]
66
AbstractFFTs = "621f4979-c628-5d54-868e-fcf4e3e8185c"
7+
BandedMatrices = "aae01518-5342-5314-be14-df237901396f"
78
FFTW = "7a1cc6ca-52ef-59f5-83cd-3a7055c09341"
89
FastGaussQuadrature = "442a2c76-b920-505d-bb47-c5924d526838"
910
FastTransforms_jll = "34b6f7d7-08f9-5794-9e10-3819e4c7e49a"
@@ -17,6 +18,7 @@ ToeplitzMatrices = "c751599d-da0a-543b-9d20-d0a503d91d24"
1718

1819
[compat]
1920
AbstractFFTs = "1.0"
21+
BandedMatrices = "1.5"
2022
FFTW = "1.7"
2123
FastGaussQuadrature = "0.4, 0.5, 1"
2224
FastTransforms_jll = "0.6.2"

src/FastTransforms.jl

+6-2
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
module FastTransforms
22

3-
using FastGaussQuadrature, FillArrays, LinearAlgebra,
3+
using BandedMatrices, FastGaussQuadrature, FillArrays, LinearAlgebra,
44
Reexport, SpecialFunctions, ToeplitzMatrices
55

66
@reexport using AbstractFFTs
@@ -19,14 +19,16 @@ import AbstractFFTs: Plan, ScaledPlan,
1919
fftshift, ifftshift, rfft_output_size, brfft_output_size,
2020
normalization
2121

22+
import BandedMatrices: bandwidths
23+
2224
import FFTW: dct, dct!, idct, idct!, plan_dct!, plan_idct!,
2325
plan_dct, plan_idct, fftwNumber
2426

2527
import FastGaussQuadrature: unweightedgausshermite
2628

2729
import FillArrays: AbstractFill, getindex_value
2830

29-
import LinearAlgebra: mul!, lmul!, ldiv!
31+
import LinearAlgebra: mul!, lmul!, ldiv!, cholesky
3032

3133
import GenericFFT: interlace # imported in downstream packages
3234

@@ -112,6 +114,8 @@ include("specialfunctions.jl")
112114
include("toeplitzplans.jl")
113115
include("toeplitzhankel.jl")
114116

117+
include("SymmetricToeplitzPlusHankel.jl")
118+
115119
# following use libfasttransforms by default
116120
for f in (:jac2jac,
117121
:lag2lag, :jac2ultra, :ultra2jac, :jac2cheb,

src/SymmetricToeplitzPlusHankel.jl

+135
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,135 @@
1+
struct SymmetricToeplitzPlusHankel{T} <: AbstractMatrix{T}
2+
v::Vector{T}
3+
n::Int
4+
end
5+
6+
function SymmetricToeplitzPlusHankel(v::Vector{T}) where T
7+
n = (length(v)+1)÷2
8+
SymmetricToeplitzPlusHankel{T}(v, n)
9+
end
10+
11+
size(A::SymmetricToeplitzPlusHankel{T}) where T = (A.n, A.n)
12+
getindex(A::SymmetricToeplitzPlusHankel{T}, i::Integer, j::Integer) where T = A.v[abs(i-j)+1] + A.v[i+j-1]
13+
14+
struct SymmetricBandedToeplitzPlusHankel{T} <: BandedMatrices.AbstractBandedMatrix{T}
15+
v::Vector{T}
16+
n::Int
17+
b::Int
18+
end
19+
20+
function SymmetricBandedToeplitzPlusHankel(v::Vector{T}, n::Integer) where T
21+
SymmetricBandedToeplitzPlusHankel{T}(v, n, length(v)-1)
22+
end
23+
24+
size(A::SymmetricBandedToeplitzPlusHankel{T}) where T = (A.n, A.n)
25+
function getindex(A::SymmetricBandedToeplitzPlusHankel{T}, i::Integer, j::Integer) where T
26+
v = A.v
27+
if abs(i-j) < length(v)
28+
if i+j-1 length(v)
29+
v[abs(i-j)+1] + v[i+j-1]
30+
else
31+
v[abs(i-j)+1]
32+
end
33+
else
34+
zero(T)
35+
end
36+
end
37+
bandwidths(A::SymmetricBandedToeplitzPlusHankel{T}) where T = (A.b, A.b)
38+
39+
#
40+
# Jac*W-W*Jac' = G*J*G'
41+
# This returns G and J, where J = [0 I; -I 0], respecting the skew-symmetry of the right-hand side.
42+
#
43+
function compute_skew_generators(A::SymmetricToeplitzPlusHankel{T}) where T
44+
v = A.v
45+
n = size(A, 1)
46+
J = [zero(T) one(T); -one(T) zero(T)]
47+
G = zeros(T, n, 2)
48+
G[n, 1] = one(T)
49+
u2 = reverse(v[2:n+1])
50+
u2[1:n-1] .+= v[n+1:2n-1]
51+
G[:, 2] .= -u2
52+
G, J
53+
end
54+
55+
function cholesky(A::SymmetricToeplitzPlusHankel{T}) where T
56+
n = size(A, 1)
57+
G, J = compute_skew_generators(A)
58+
L = zeros(T, n, n)
59+
r = A[:, 1]
60+
r2 = zeros(T, n)
61+
l = zeros(T, n)
62+
v = zeros(T, n)
63+
col1 = zeros(T, n)
64+
STPHcholesky!(L, G, r, r2, l, v, col1, n)
65+
return Cholesky(L, 'L', 0)
66+
end
67+
68+
function STPHcholesky!(L::Matrix{T}, G, r, r2, l, v, col1, n) where T
69+
@inbounds @simd for k in 1:n-1
70+
x = sqrt(r[1])
71+
for j in 1:n-k+1
72+
L[j+k-1, k] = l[j] = r[j]/x
73+
end
74+
for j in 1:n-k+1
75+
v[j] = G[j, 1]*G[1,2]-G[j,2]*G[1,1]
76+
end
77+
F12 = k == 1 ? T(2) : T(1)
78+
r2[1] = (r[2] - v[1])/F12
79+
for j in 2:n-k
80+
r2[j] = (r[j+1]+r[j-1] + r[1]*col1[j] - col1[1]*r[j] - v[j])/F12
81+
end
82+
r2[n-k+1] = (r[n-k] + r[1]*col1[n-k+1] - col1[1]*r[n-k+1] - v[n-k+1])/F12
83+
cst = r[2]/x
84+
for j in 1:n-k
85+
r[j] = r2[j+1] - cst*l[j+1]
86+
end
87+
for j in 1:n-k
88+
col1[j] = -F12/x*l[j+1]
89+
end
90+
c1 = G[1, 1]
91+
c2 = G[1, 2]
92+
for j in 1:n-k
93+
G[j, 1] = G[j+1, 1] - l[j+1]*c1/x
94+
G[j, 2] = G[j+1, 2] - l[j+1]*c2/x
95+
end
96+
end
97+
L[n, n] = sqrt(r[1])
98+
end
99+
100+
function cholesky(A::SymmetricBandedToeplitzPlusHankel{T}) where T
101+
n = size(A, 1)
102+
b = A.b
103+
R = BandedMatrix{T}(undef, (n, n), (0, bandwidth(A, 2)))
104+
r = A[1:b+2, 1]
105+
r2 = zeros(T, b+3)
106+
l = zeros(T, b+3)
107+
col1 = zeros(T, b+2)
108+
SBTPHcholesky!(R, r, r2, l, col1, n, b)
109+
return Cholesky(R, 'U', 0)
110+
end
111+
112+
function SBTPHcholesky!(R::BandedMatrix{T}, r, r2, l, col1, n, b) where T
113+
@inbounds @simd for k in 1:n
114+
x = sqrt(r[1])
115+
for j in 1:b+1
116+
l[j] = r[j]/x
117+
end
118+
for j in 1:min(n-k+1, b+1)
119+
R[k, j+k-1] = l[j]
120+
end
121+
F12 = k == 1 ? T(2) : T(1)
122+
r2[1] = r[2]/F12
123+
for j in 2:b+1
124+
r2[j] = (r[j+1]+r[j-1] + r[1]*col1[j] - col1[1]*r[j])/F12
125+
end
126+
r2[b+2] = (r[b+1] + r[1]*col1[b+2] - col1[1]*r[b+2])/F12
127+
cst = r[2]/x
128+
for j in 1:b+2
129+
r[j] = r2[j+1] - cst*l[j+1]
130+
end
131+
for j in 1:b+2
132+
col1[j] = -F12/x*l[j+1]
133+
end
134+
end
135+
end

test/runtests.jl

+2-1
Original file line numberDiff line numberDiff line change
@@ -10,4 +10,5 @@ include("gaunttests.jl")
1010
include("hermitetests.jl")
1111
include("clenshawtests.jl")
1212
include("toeplitzplanstests.jl")
13-
include("toeplitzhankeltests.jl")
13+
include("toeplitzhankeltests.jl")
14+
include("symmetrictoeplitzplushankeltests.jl")
+39
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,39 @@
1+
using BandedMatrices, FastTransforms, LinearAlgebra, ToeplitzMatrices, Test
2+
3+
import FastTransforms: SymmetricToeplitzPlusHankel, SymmetricBandedToeplitzPlusHankel
4+
5+
@testset "SymmetricToeplitzPlusHankel" begin
6+
n = 128
7+
for T in (Float32, Float64, BigFloat)
8+
μ = -FastTransforms.chebyshevlogmoments1(T, 2n-1)
9+
μ[1] += 1
10+
W = SymmetricToeplitzPlusHankel/2)
11+
SMW = Symmetric(Matrix(W))
12+
@test W SymmetricToeplitz(μ[1:(length(μ)+1)÷2]/2) + Hankel/2)
13+
L = cholesky(W).L
14+
R = cholesky(SMW).U
15+
@test L*L' W
16+
@test L' R
17+
end
18+
end
19+
20+
@testset "SymmetricBandedToeplitzPlusHankel" begin
21+
n = 1024
22+
for T in (Float32, Float64)
23+
μ = T[1.875; 0.00390625; 0.5; 0.0009765625; 0.0625]
24+
W = SymmetricBandedToeplitzPlusHankel/2, n)
25+
SBW = Symmetric(BandedMatrix(W))
26+
W1 = SymmetricToeplitzPlusHankel([μ/2; zeros(2n-1-length(μ))])
27+
SMW = Symmetric(Matrix(W))
28+
U = cholesky(SMW).U
29+
L = cholesky(W1).L
30+
UB = cholesky(SBW).U
31+
R = cholesky(W).U
32+
@test L*L' W
33+
@test UB'UB W
34+
@test R'R W
35+
@test UB U
36+
@test L' U
37+
@test R U
38+
end
39+
end

0 commit comments

Comments
 (0)