Skip to content

Commit 14507a0

Browse files
authored
Support vec with RectPolynomial and fix sum ordering (#208)
* Support vec with RectPolynomial and fix sum ordering * reshape and vec * Update test_rect.jl
1 parent 22234cc commit 14507a0

File tree

4 files changed

+62
-14
lines changed

4 files changed

+62
-14
lines changed

Project.toml

Lines changed: 2 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -30,7 +30,7 @@ BandedMatrices = "1"
3030
BlockArrays = "1.0"
3131
BlockBandedMatrices = "0.13"
3232
ClassicalOrthogonalPolynomials = "0.15"
33-
ContinuumArrays = "0.19"
33+
ContinuumArrays = "0.19.5"
3434
DomainSets = "0.7"
3535
FastTransforms = "0.17"
3636
FillArrays = "1.0"
@@ -39,7 +39,7 @@ InfiniteArrays = "0.15"
3939
InfiniteLinearAlgebra = "0.9, 0.10"
4040
LazyArrays = "2.3.1"
4141
LazyBandedMatrices = "0.11.3"
42-
QuasiArrays = "0.12"
42+
QuasiArrays = "0.12.2"
4343
RecurrenceRelationships = "0.2"
4444
SpecialFunctions = "1, 2"
4545
StaticArrays = "1"

src/MultivariateOrthogonalPolynomials.jl

Lines changed: 3 additions & 3 deletions
Original file line numberDiff line numberDiff line change
@@ -10,8 +10,8 @@ import Base: axes, in, ==, +, -, /, *, ^, \, copy, copyto!, OneTo, getindex, siz
1010
import Base.Broadcast: Broadcasted, broadcasted, DefaultArrayStyle
1111
import DomainSets: boundary, EuclideanDomain
1212

13-
import QuasiArrays: LazyQuasiMatrix, LazyQuasiArrayStyle, domain
14-
import ContinuumArrays: @simplify, Weight, weight, grid, plotgrid, TransformFactorization, ExpansionLayout, plotvalues, unweighted, plan_transform, checkpoints, transform_ldiv, AbstractBasisLayout, basis_axes, Inclusion, grammatrix, weaklaplacian, layout_broadcasted, laplacian, abslaplacian, laplacian_axis, abslaplacian_axis, diff_layout, operatororder, broadcastbasis
13+
import QuasiArrays: LazyQuasiMatrix, LazyQuasiArrayStyle, domain, vec_layout, reshape_layout
14+
import ContinuumArrays: @simplify, Weight, weight, grid, plotgrid, TransformFactorization, ExpansionLayout, plotvalues, unweighted, plan_transform, checkpoints, transform_ldiv, AbstractBasisLayout, basis_axes, Inclusion, grammatrix, weaklaplacian, layout_broadcasted, laplacian, abslaplacian, laplacian_axis, abslaplacian_axis, diff_layout, operatororder, broadcastbasis, KronExpansionLayout
1515

1616
import ArrayLayouts: MemoryLayout, sublayout, sub_materialize
1717
import BlockArrays: block, blockindex, BlockSlice, viewblock, blockcolsupport, AbstractBlockStyle, BlockStyle
@@ -21,7 +21,7 @@ import LazyArrays: arguments, paddeddata, LazyArrayStyle, LazyLayout, PaddedLayo
2121
import LazyBandedMatrices: LazyBandedBlockBandedLayout, AbstractBandedBlockBandedLayout, AbstractLazyBandedBlockBandedLayout, _krontrav_axes, DiagTravLayout, invdiagtrav, ApplyBandedBlockBandedLayout, krontrav
2222
import InfiniteArrays: InfiniteCardinal, OneToInf
2323

24-
import ClassicalOrthogonalPolynomials: jacobimatrix, Weighted, orthogonalityweight, HalfWeighted, WeightedBasis, pad, recurrencecoefficients, clenshaw, weightedgrammatrix, Clenshaw
24+
import ClassicalOrthogonalPolynomials: jacobimatrix, Weighted, orthogonalityweight, HalfWeighted, WeightedBasis, pad, recurrencecoefficients, clenshaw, weightedgrammatrix, Clenshaw, OPLayout
2525
import HarmonicOrthogonalPolynomials: BivariateOrthogonalPolynomial, MultivariateOrthogonalPolynomial, Plan,
2626
AngularMomentum, angularmomentum, BlockOneTo, BlockRange1, interlace,
2727
MultivariateOPLayout, AbstractMultivariateOPLayout, MAX_PLOT_BLOCKS

src/rect.jl

Lines changed: 21 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -27,6 +27,8 @@ const RectPolynomial{T, PP} = KronPolynomial{2, T, PP}
2727

2828

2929
axes(P::KronPolynomial) = (Inclusion(×(map(domain, axes.(P.args, 1))...)), _krontrav_axes(axes.(P.args, 2)...))
30+
31+
3032
function getindex(P::RectPolynomial{T}, xy::StaticVector{2}, Jj::BlockIndex{1})::T where T
3133
a,b = P.args
3234
J,j = Int(block(Jj)),blockindex(Jj)
@@ -166,7 +168,7 @@ end
166168

167169
function Base._sum(P::RectPolynomial, dims::Int)
168170
@assert dims == 1
169-
KronTrav(sum.(P.args; dims=1)...)
171+
KronTrav(reverse(sum.(P.args; dims=1))...)
170172
end
171173

172174
## multiplication
@@ -183,4 +185,21 @@ function layout_broadcasted(::Tuple{ExpansionLayout{KronOPLayout{2}},KronOPLayou
183185
end
184186

185187

186-
broadcastbasis(::typeof(+), A::KronPolynomial, B::KronPolynomial) = KronPolynomial(broadcastbasis.(+, A.args, B.args)...)
188+
broadcastbasis(::typeof(+), A::KronPolynomial, B::KronPolynomial) = KronPolynomial(broadcastbasis.(+, A.args, B.args)...)
189+
190+
191+
192+
####
193+
# reshape/vec
194+
####
195+
196+
function reshape_layout(::ExpansionLayout{KronOPLayout{2}}, f)
197+
A,B = basis(f).args
198+
A*invdiagtrav(coefficients(f))*B'
199+
end
200+
201+
# TODO: generalise
202+
function vec_layout(::KronExpansionLayout{OPLayout, OPLayout}, f)
203+
A,C,Bc = arguments(f)
204+
RectPolynomial(A, Bc') * DiagTrav(C)
205+
end

test/test_rect.jl

Lines changed: 36 additions & 7 deletions
Original file line numberDiff line numberDiff line change
@@ -1,7 +1,7 @@
1-
using MultivariateOrthogonalPolynomials, ClassicalOrthogonalPolynomials, StaticArrays, LinearAlgebra, BlockArrays, FillArrays, Base64, LazyBandedMatrices, Test
1+
using MultivariateOrthogonalPolynomials, ClassicalOrthogonalPolynomials, StaticArrays, LinearAlgebra, BlockArrays, FillArrays, Base64, LazyBandedMatrices, ArrayLayouts, Test
22
using ClassicalOrthogonalPolynomials: expand, coefficients, recurrencecoefficients
33
using MultivariateOrthogonalPolynomials: weaklaplacian, ClenshawKron
4-
using ContinuumArrays: plotgridvalues
4+
using ContinuumArrays: plotgridvalues, ExpansionLayout
55
using Base: oneto
66

77
@testset "RectPolynomial" begin
@@ -148,11 +148,12 @@ using Base: oneto
148148
end
149149

150150
@testset "sum" begin
151-
P = RectPolynomial(Legendre(),Legendre())
152-
p₀ = expand(P, 𝐱 -> 1)
153-
@test sum(p₀) 4.0
154-
f = expand(P, splat((x,y) -> exp(cos(x^2*y))))
155-
@test sum(f) 10.546408460894801 # empirical
151+
for P in (RectPolynomial(Legendre(),Legendre()), RectPolynomial(Legendre(),Chebyshev()))
152+
p₀ = expand(P, 𝐱 -> 1)
153+
@test sum(p₀) 4.0
154+
f = expand(P, splat((x,y) -> exp(cos(x^2*y))))
155+
@test sum(f) 10.546408460894801 # empirical
156+
end
156157
end
157158

158159
@testset "diff" begin
@@ -219,4 +220,32 @@ using Base: oneto
219220

220221
𝐜 = expand(RectPolynomial(Legendre(),Jacobi(1,0)),splat((x,y) -> cos(x*sin(y))))
221222
end
223+
224+
@testset "reshape/vec" begin
225+
P = RectPolynomial(Legendre(),Chebyshev())
226+
f = expand(P, splat((x,y) -> cos((x-0.1)*exp(y))))
227+
F = reshape(f)
228+
@test F[0.1,0.2] f[SVector(0.1,0.2)]
229+
@test vec(F)[SVector(0.1,0.2)] f[SVector(0.1,0.2)]
230+
231+
g = F[:,0.2]
232+
h = F[0.1,:]
233+
@test MemoryLayout(g) isa ExpansionLayout
234+
@test MemoryLayout(h) isa ExpansionLayout
235+
@test g[0.1] f[SVector(0.1,0.2)]
236+
@test h[0.2] f[SVector(0.1,0.2)]
237+
238+
@test sum(F; dims=1)[1,0.2] exp(-0.2)*(sin(0.9exp(0.2)) + sin(1.1exp(0.2)))
239+
# TODO: should be matrix but isn't because of InfiniteArrays/src/reshapedarray.jl:77
240+
@test_broken sum(F; dims=2)[0.1,1] 2
241+
@test sum(F; dims=2)[0.1] 2
242+
end
243+
244+
@testset "sample" begin
245+
P = RectPolynomial(Legendre(),Legendre())
246+
f = expand(P, splat((x,y) -> exp(x*cos(y-0.1))))
247+
F = reshape(f)
248+
@test sum(F; dims=1)[1,0.2] 2.346737615950585
249+
@test sum(F; dims=2)[0.1,1] 2.1748993079723618
250+
end
222251
end

0 commit comments

Comments
 (0)