@@ -7,13 +7,14 @@ import Base: getindex, axes, size, \, /, *, +, -, summary, ==, copy, sum, unsafe
7
7
8
8
import ArrayLayouts: MemoryLayout, ldiv, diagonaldata, subdiagonaldata, supdiagonaldata
9
9
import BandedMatrices: bandwidths, AbstractBandedMatrix, BandedLayout, _BandedMatrix
10
- import LazyArrays: resizedata!, paddeddata, CachedVector, CachedMatrix, CachedAbstractVector, LazyMatrix, LazyVector, arguments, ApplyLayout, colsupport, AbstractCachedVector,
10
+ import LazyArrays: resizedata!, paddeddata, CachedVector, CachedMatrix, CachedAbstractVector, LazyMatrix, LazyVector, arguments, ApplyLayout, colsupport, AbstractCachedVector, ApplyArray,
11
11
AccumulateAbstractVector, LazyVector, AbstractCachedMatrix, BroadcastLayout
12
12
import ClassicalOrthogonalPolynomials: OrthogonalPolynomial, recurrencecoefficients, jacobimatrix, normalize, _p0, UnitInterval, orthogonalityweight, NormalizedOPLayout,
13
13
Bidiagonal, Tridiagonal, SymTridiagonal, symtridiagonalize, normalizationconstant, LanczosPolynomial,
14
14
OrthogonalPolynomialRatio, Weighted, WeightLayout, UnionDomain, oneto, WeightedBasis, HalfWeighted,
15
15
golubwelsch, AbstractOPLayout, weight, WeightedOPLayout
16
16
import SingularIntegrals: Associated, associated, Hilbert
17
+
17
18
import InfiniteArrays: OneToInf, InfUnitRange
18
19
import ContinuumArrays: basis, Weight, @simplify , AbstractBasisLayout, BasisLayout, MappedBasisLayout, grid, plotgrid, _equals, ExpansionLayout
19
20
import FillArrays: SquareEye
@@ -94,20 +95,13 @@ RaisedOP(Q, y::Number) = RaisedOP(Q, OrthogonalPolynomialRatio(Q,y))
94
95
function jacobimatrix (P:: RaisedOP{T} ) where T
95
96
ℓ = P. ℓ
96
97
X = jacobimatrix (P. Q)
97
- a,b = diagonaldata (X), supdiagonaldata (X)
98
+ a,b = view (X, band ( 0 )), view (X, band ( 1 )) # a,b = diagonaldata(X), supdiagonaldata(X)
98
99
# non-normalized lower diag of Jacobi
99
100
v = Vcat (zero (T),b .* ℓ)
100
101
c = BroadcastVector ((ℓ,a,b,sa,v) -> ℓ* a + b - b* ℓ^ 2 - sa* ℓ + ℓ* v, ℓ, a, b, a[2 : ∞], v)
101
102
Tridiagonal (c, BroadcastVector ((ℓ,a,b,v) -> a - b * ℓ + v, ℓ,a,b,v), b)
102
103
end
103
104
104
-
105
-
106
-
107
-
108
-
109
-
110
-
111
105
"""
112
106
the bands of the Jacobi matrix
113
107
"""
@@ -142,7 +136,6 @@ SemiclassicalJacobiBand{dv}(a,b,ℓ) where dv = SemiclassicalJacobiBand{dv,promo
142
136
copy (r:: SemiclassicalJacobiBand ) = r # immutable
143
137
144
138
145
-
146
139
"""
147
140
SemiclassicalJacobi(t, a, b, c)
148
141
@@ -177,7 +170,7 @@ function semiclassical_jacobimatrix(t, a, b, c)
177
170
if a == 0 && b == 0 && c == - 1
178
171
# for this special case we can generate the Jacobi operator explicitly
179
172
N = (1 : ∞)
180
- α = T .( neg1c_αcfs (t) ) # cached coefficients
173
+ α = neg1c_αcfs (t) # cached coefficients
181
174
A = Vcat ((α[1 ]+ 1 )/ 2 , - N./ (N.* 4 .- 2 ). * α .+ (N.+ 1 ). / (N.* 4 .+ 2 ). * α[2 : end ]. + 1 / 2 )
182
175
C = - (N). / (N.* 4 .- 2 )
183
176
B = Vcat ((α[1 ]^ 2 * 3 - α[1 ]* α[2 ]* 2 - 1 )/ 6 , - (N). / (N.* 4 .+ 2 ). * α[2 : end ]. / α)
@@ -187,10 +180,15 @@ function semiclassical_jacobimatrix(t, a, b, c)
187
180
P = jacobi (b, a, UnitInterval {T} ())
188
181
iszero (c) && return symtridiagonalize (jacobimatrix (P))
189
182
x = axes (P,1 )
190
- jacobimatrix (LanczosPolynomial (@. (x^ a * (1 - x)^ b * (t- x)^ c), P))
183
+ if iseven (c) # QR case
184
+ return qr_jacobimatrix (x-> (t.- x). ^ (c÷ 2 ), Normalized (jacobi (a,b,UnitInterval {T} ())))
185
+ elseif isodd (c)
186
+ return semiclassical_jacobimatrix (SemiclassicalJacobi (t, a, b, c- 1 ), a, b, c)
187
+ else # if c is not an even integer, use Lanczos for now
188
+ return jacobimatrix (LanczosPolynomial (@. (x^ a * (1 - x)^ b * (t- x)^ c), P))
189
+ end
191
190
end
192
191
193
-
194
192
function symraised_jacobimatrix (Q, y)
195
193
ℓ = OrthogonalPolynomialRatio (Q,y)
196
194
X = jacobimatrix (Q)
@@ -199,9 +197,35 @@ function symraised_jacobimatrix(Q, y)
199
197
end
200
198
201
199
function semiclassical_jacobimatrix (Q:: SemiclassicalJacobi , a, b, c)
202
- if iszero (a) && iszero (b) && c == - 1 # (a,b,c) = (0,0,-1) special case
200
+ # special cases
201
+ if iszero (a) && iszero (b) && c == - 1 # (a,b,c) = (0,0,-1) special case
203
202
semiclassical_jacobimatrix (Q. t, a, b, c)
204
- elseif a == Q. a+ 1 && b == Q. b && c == Q. c # raising by 1
203
+ elseif iszero (c) # classical Jacobi polynomial special case
204
+ jacobimatrix (SemiclassicalJacobi (Q. t, Q. a, Q. b, Q. c))
205
+ elseif a == Q. a && b == Q. b && c == Q. c # same basis
206
+ jacobimatrix (Q)
207
+ end
208
+
209
+ # QR decomposition approach for raising
210
+ Δa = a- Q. a
211
+ Δb = b- Q. b
212
+ Δc = c- Q. c
213
+ if iseven (a- Q. a) && iseven (b- Q. b) && iseven (c- Q. c) && Δa≥ 0 && Δb≥ 0 && Δc≥ 0
214
+ M = Diagonal (Ones (∞))
215
+ if ! iszero (Δa)
216
+ M = ApplyArray (* ,M,(Q. X)^ ((a- Q. a)÷ 2 ))
217
+ end
218
+ if ! iszero (Δb)
219
+ M = ApplyArray (* ,M,(I- Q. X)^ ((b- Q. b)÷ 2 ))
220
+ end
221
+ if ! iszero (Δc)
222
+ M = ApplyArray (* ,M,(Q. t* I- Q. X)^ ((c- Q. c)÷ 2 ))
223
+ end
224
+ return qr_jacobimatrix (M,Q,false )
225
+ end
226
+
227
+ # Fallback from decomposition methods
228
+ if a == Q. a+ 1 && b == Q. b && c == Q. c # raising by 1
205
229
symraised_jacobimatrix (Q, 0 )
206
230
elseif a == Q. a && b == Q. b+ 1 && c == Q. c
207
231
symraised_jacobimatrix (Q, 1 )
@@ -213,20 +237,18 @@ function semiclassical_jacobimatrix(Q::SemiclassicalJacobi, a, b, c)
213
237
symlowered_jacobimatrix (Q, :a )
214
238
elseif a == Q. a && b == Q. b- 1 && c == Q. c
215
239
symlowered_jacobimatrix (Q, :b )
216
- elseif a > Q. a # iterative raising
240
+ elseif a > Q. a # iterative raising by 1
217
241
semiclassical_jacobimatrix (SemiclassicalJacobi (Q. t, Q. a+ 1 , Q. b, Q. c, Q), a, b, c)
218
242
elseif b > Q. b
219
243
semiclassical_jacobimatrix (SemiclassicalJacobi (Q. t, Q. a, Q. b+ 1 , Q. c, Q), a, b, c)
220
244
elseif c > Q. c
221
245
semiclassical_jacobimatrix (SemiclassicalJacobi (Q. t, Q. a, Q. b, Q. c+ 1 , Q), a, b, c)
222
- elseif b < Q. b # iterative lowering
246
+ elseif b < Q. b # iterative lowering by 1
223
247
semiclassical_jacobimatrix (SemiclassicalJacobi (Q. t, Q. a, Q. b- 1 , Q. c, Q), a, b, c)
224
248
elseif c < Q. c
225
249
semiclassical_jacobimatrix (SemiclassicalJacobi (Q. t, Q. a, Q. b, Q. c- 1 , Q), a, b, c)
226
250
elseif a < Q. a
227
251
semiclassical_jacobimatrix (SemiclassicalJacobi (Q. t, Q. a- 1 , Q. b, Q. c, Q), a, b, c)
228
- elseif a == Q. a && b == Q. b && c == Q. c # same basis
229
- jacobimatrix (Q)
230
252
else
231
253
error (" Not Implemented" )
232
254
end
@@ -292,6 +314,46 @@ function semijacobi_ldiv(P::SemiclassicalJacobi, Q)
292
314
(P \ R) * _p0 (R̃) * (R̃ \ Q)
293
315
end
294
316
317
+ # returns conversion operator from SemiclassicalJacobi P to SemiclassicalJacobi Q.
318
+ function semijacobi_ldiv (Q:: SemiclassicalJacobi ,P:: SemiclassicalJacobi )
319
+ @assert Q. t ≈ P. t
320
+ # For polynomial modifications, use Cholesky. Use Lanzos otherwise.
321
+ M = Diagonal (Ones (∞))
322
+ (Q == P) && return M
323
+ Δa = Q. a- P. a
324
+ Δb = Q. b- P. b
325
+ Δc = Q. c- P. c
326
+ if iseven (Δa) && iseven (Δb) && iseven (Δc)
327
+ if ! iszero (Q. a- P. a)
328
+ M = ApplyArray (* ,M,(P. X)^ ((Q. a- P. a)÷ 2 ))
329
+ end
330
+ if ! iszero (Q. b- P. b)
331
+ M = ApplyArray (* ,M,(I- P. X)^ ((Q. b- P. b)÷ 2 ))
332
+ end
333
+ if ! iszero (Q. c- P. c)
334
+ M = ApplyArray (* ,M,(Q. t* I- P. X)^ ((Q. c- P. c)÷ 2 ))
335
+ end
336
+ K = qr (M). R
337
+ return ApplyArray (* , Diagonal (sign .(view (K,band (0 ))).* Fill (abs .(1 / K[1 ]),∞)), K) # match our normalization choice P_0(x) = 1
338
+ elseif isinteger (Δa) && isinteger (Δb) && isinteger (Δc)
339
+ if ! iszero (Q. a- P. a)
340
+ M = ApplyArray (* ,M,(P. X)^ (Q. a- P. a))
341
+ end
342
+ if ! iszero (Q. b- P. b)
343
+ M = ApplyArray (* ,M,(I- P. X)^ (Q. b- P. b))
344
+ end
345
+ if ! iszero (Q. c- P. c)
346
+ M = ApplyArray (* ,M,(Q. t* I- P. X)^ (Q. c- P. c))
347
+ end
348
+ K = cholesky (Symmetric (M)). U
349
+ return ApplyArray (* , Diagonal (Fill (1 / K[1 ],∞)), K) # match our normalization choice P_0(x) = 1
350
+ else # fallback option for non-integer weight modification
351
+ R = SemiclassicalJacobi (P. t, mod (P. a,- 1 ), mod (P. b,- 1 ), mod (P. c,- 1 ))
352
+ R̃ = toclassical (R)
353
+ return (P \ R) * _p0 (R̃) * (R̃ \ Q)
354
+ end
355
+ end
356
+
295
357
struct SemiclassicalJacobiLayout <: AbstractOPLayout end
296
358
MemoryLayout (:: Type{<:SemiclassicalJacobi} ) = SemiclassicalJacobiLayout ()
297
359
@@ -323,6 +385,7 @@ function copy(L::Ldiv{SemiclassicalJacobiLayout,SemiclassicalJacobiLayout})
323
385
(inv (M_Q) * L' ) * M_P
324
386
end
325
387
388
+ \ (A:: SemiclassicalJacobi , B:: SemiclassicalJacobi ) = semijacobi_ldiv (A, B)
326
389
\ (A:: LanczosPolynomial , B:: SemiclassicalJacobi ) = semijacobi_ldiv (A, B)
327
390
\ (A:: SemiclassicalJacobi , B:: LanczosPolynomial ) = semijacobi_ldiv (A, B)
328
391
function \ (w_A:: WeightedSemiclassicalJacobi , w_B:: WeightedSemiclassicalJacobi )
@@ -456,7 +519,27 @@ function SemiclassicalJacobiFamily{T}(data::Vector, t, a, b, c) where T
456
519
end
457
520
458
521
SemiclassicalJacobiFamily (t, a, b, c) = SemiclassicalJacobiFamily {float(promote_type(typeof(t),eltype(a),eltype(b),eltype(c)))} (t, a, b, c)
459
- SemiclassicalJacobiFamily {T} (t, a, b, c) where T = SemiclassicalJacobiFamily {T} ([SemiclassicalJacobi {T} (t, first (a), first (b), first (c))], t, a, b, c)
522
+ function SemiclassicalJacobiFamily {T} (t, a, b, c) where T
523
+ ℓa = length (a)
524
+ ℓb = length (b)
525
+ ℓc = length (c)
526
+ # We need to start with a hierarchy containing two entries
527
+ if ℓa > 1 && ℓb > 1 && ℓc > 1
528
+ return SemiclassicalJacobiFamily {T} ([SemiclassicalJacobi {T} (t, first (a), first (b), first (c)),SemiclassicalJacobi {T} (t, a[2 ], b[2 ], c[2 ])], t, a, b, c)
529
+ elseif ℓb > 1 && ℓc > 1
530
+ return SemiclassicalJacobiFamily {T} ([SemiclassicalJacobi {T} (t, first (a), first (b), first (c)),SemiclassicalJacobi {T} (t, first (a), b[2 ], c[2 ])], t, a, b, c)
531
+ elseif ℓa > 1 && ℓc > 1
532
+ return SemiclassicalJacobiFamily {T} ([SemiclassicalJacobi {T} (t, first (a), first (b), first (c)),SemiclassicalJacobi {T} (t, a[2 ], first (b), c[2 ])], t, a, b, c)
533
+ elseif ℓa > 1 && ℓb > 1
534
+ return SemiclassicalJacobiFamily {T} ([SemiclassicalJacobi {T} (t, first (a), first (b), first (c)),SemiclassicalJacobi {T} (t, a[2 ], b[2 ], first (c))], t, a, b, c)
535
+ elseif ℓc > 1
536
+ return SemiclassicalJacobiFamily {T} ([SemiclassicalJacobi {T} (t, first (a), first (b), first (c)),SemiclassicalJacobi {T} (t, first (a), first (b), c[2 ])], t, a, b, c)
537
+ elseif ℓa > 1
538
+ return SemiclassicalJacobiFamily {T} ([SemiclassicalJacobi {T} (t, first (a), first (b), first (c)),SemiclassicalJacobi {T} (t, a[2 ], first (b), first (c))], t, a, b, c)
539
+ else # if ℓb > 1
540
+ return SemiclassicalJacobiFamily {T} ([SemiclassicalJacobi {T} (t, first (a), first (b), first (c)),SemiclassicalJacobi {T} (t, first (a), b[2 ], first (c))], t, a, b, c)
541
+ end
542
+ end
460
543
461
544
Base. broadcasted (:: Type{SemiclassicalJacobi} , t:: Number , a:: Number , b:: Number , c:: Number ) = SemiclassicalJacobi (t, a, b, c)
462
545
Base. broadcasted (:: Type{SemiclassicalJacobi{T}} , t:: Number , a:: Number , b:: Number , c:: Number ) where T = SemiclassicalJacobi {T} (t, a, b, c)
@@ -472,7 +555,7 @@ _broadcast_getindex(a::Number,k) = a
472
555
function LazyArrays. cache_filldata! (P:: SemiclassicalJacobiFamily , inds:: AbstractUnitRange )
473
556
t,a,b,c = P. t,P. a,P. b,P. c
474
557
for k in inds
475
- P. data[k] = SemiclassicalJacobi (t, _broadcast_getindex (a,k), _broadcast_getindex (b,k), _broadcast_getindex (c,k), P. data[k- 1 ])
558
+ P. data[k] = SemiclassicalJacobi (t, _broadcast_getindex (a,k), _broadcast_getindex (b,k), _broadcast_getindex (c,k), P. data[k- 2 ])
476
559
end
477
560
P
478
561
end
0 commit comments