Skip to content

Commit 36b0dbe

Browse files
Merge pull request #466 from cgarling/alloc-checks
Test allocations
2 parents 8e6fbe4 + b119ee6 commit 36b0dbe

3 files changed

Lines changed: 97 additions & 8 deletions

File tree

src/interpolation_utils.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -38,7 +38,7 @@ function spline_coefficients!(N, d, k, u::Number)
3838
N[end] = one(u)
3939
return length(N):length(N)
4040
else
41-
i = findfirst(x -> x > u, k) - 1
41+
i = findfirst(x -> x > u, k)::Int - 1
4242
N[i] = one(u)
4343
for deg in 1:d
4444
N[i - deg] = (k[i + 1] - u) / (k[i + 1] - k[i - deg + 1]) * N[i - deg + 1]

src/parameter_caches.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -61,7 +61,7 @@ function smoothed_constant_interpolation_parameters(
6161
else
6262
d = isone(idx) ? min(t[2] - t[1], 2d_max) / 2 :
6363
min(t[end] - t[end - 1], 2d_max) / 2
64-
d, zero(one(eltype(u)) / 2)
64+
d, zero(first(u) / 2)
6565
end
6666
else
6767
min(t[idx] - t[idx - 1], t[idx + 1] - t[idx], 2d_max) / 2, (u[idx] - u[idx - 1]) / 2

test/interpolation_tests.jl

Lines changed: 95 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -90,7 +90,7 @@ end
9090
@test @inferred(LinearInterpolation(
9191
u_s, t; extrapolation = ExtrapolationType.Extension)) isa LinearInterpolation
9292
A_s = LinearInterpolation(u_s, t; extrapolation = ExtrapolationType.Extension)
93-
for _x in (0, 5.5, 11)
93+
for x in (0, 5.5, 11)
9494
@test A(x) == A_s(x)
9595
end
9696
@test A_s(0) isa SVector{length(y)}
@@ -331,6 +331,7 @@ end
331331
@test @inferred(output_dim(A)) == 1
332332
@test @inferred(output_size(A)) == (2,)
333333

334+
# Vector{Vector} interpolation test
334335
u_ = [1.0, 4.0, 9.0, 16.0]' .* ones(5)
335336
u = [u_[:, i] for i in 1:size(u_, 2)]
336337
A = @inferred(QuadraticInterpolation(u, t; extrapolation = ExtrapolationType.Extension))
@@ -341,7 +342,18 @@ end
341342
@test A(5.0) == 25.0 * ones(5)
342343
@test @inferred(output_dim(A)) == 1
343344
@test @inferred(output_size(A)) == (5,)
345+
# Test allocation-free interpolation with Vector{StaticArrays.SVector}
346+
u_s = [convert(SVector{length(u[1])}, i) for i in u]
347+
@test @inferred(QuadraticInterpolation(
348+
u_s, t; extrapolation = ExtrapolationType.Extension)) isa QuadraticInterpolation
349+
A_s = QuadraticInterpolation(u_s, t; extrapolation = ExtrapolationType.Extension)
350+
for x in (0, 1.5, 2.5, 3.5, 5.0)
351+
@test A(x) == A_s(x)
352+
end
353+
@test A_s(0) isa SVector{length(u[1])}
354+
@test_nowarn test_allocs(A_s, 0)
344355

356+
# Vector{Matrix} interpolation test
345357
u = [repeat(u[i], 1, 3) for i in 1:4]
346358
@test_broken @inferred(QuadraticInterpolation(
347359
u, t; extrapolation = ExtrapolationType.Extension)) isa QuadraticInterpolation
@@ -568,6 +580,16 @@ end
568580
test_cached_index(A)
569581
@test @inferred(output_dim(A)) == 1
570582
@test @inferred(output_size(A)) == (2,)
583+
584+
# Test allocation-free interpolation with StaticArrays
585+
u_s = [convert(SVector{length(first(u))}, i) for i in u]
586+
A_s = @inferred(ConstantInterpolation(
587+
u_s, t; extrapolation = ExtrapolationType.Extension))
588+
for x in 0.5:0.5:4.5
589+
@test A(x) == A_s(x)
590+
end
591+
@test A_s(0) isa SVector{length(first(u))}
592+
@test_nowarn test_allocs(A_s, 0)
571593
end
572594

573595
@testset "Vector of Matrices case" for u in [
@@ -647,6 +669,30 @@ end
647669
@test A(1.9) == u[1]
648670
@test A(3.1) == u[2]
649671
@test A(2.5) (u[1] + u[2]) / 2
672+
673+
u_ = u' .* ones(5) # [u for i in 1:length(t)]
674+
uv = [u_[:, i] for i in 1:size(u_, 2)]
675+
# Test Vector{Vector} interpolation
676+
A = SmoothedConstantInterpolation(uv, t; d_max)
677+
@test A(1.9) == uv[1]
678+
@test A(3.1) == uv[2]
679+
@test A(2.5) (uv[1] + uv[2]) / 2
680+
# Test allocation-free interpolation with Vector{StaticArrays.SVector}
681+
u_s = [convert(SVector{length(uv[1])}, i) for i in uv]
682+
@test_broken @inferred(SmoothedConstantInterpolation(
683+
u_s, t; d_max)) isa SmoothedConstantInterpolation
684+
A_s = SmoothedConstantInterpolation(u_s, t; d_max)
685+
for x in (1.9, 3.1, 2.5)
686+
@test A(x) == A_s(x)
687+
end
688+
@test A_s(1.9) isa SVector{length(uv[1])}
689+
@test_nowarn test_allocs(A_s, 1.9)
690+
# Test Vector{Matrix} interpolation
691+
um = [repeat(uv[i], 1, 3) for i in 1:length(t)]
692+
A = SmoothedConstantInterpolation(um, t; d_max)
693+
@test A(1.9) == u[1] * ones(5, 3)
694+
@test A(3.1) == u[2] * ones(5, 3)
695+
@test A(2.5) ((u[1] + u[2]) / 2) * ones(5, 3)
650696
end
651697

652698
@testset "QuadraticSpline Interpolation" begin
@@ -671,6 +717,7 @@ end
671717
@test @inferred(output_dim(A)) == 0
672718
@test @inferred(output_size(A)) == ()
673719

720+
# Test Vector{Vector} interpolation
674721
u_ = [0.0, 1.0, 3.0]' .* ones(4)
675722
u = [u_[:, i] for i in 1:size(u_, 2)]
676723
A = QuadraticSpline(u, t; extrapolation = ExtrapolationType.Extension)
@@ -680,7 +727,16 @@ end
680727
@test A(2.0) == P₁(2.0) * ones(4)
681728
@test @inferred(output_dim(A)) == 1
682729
@test @inferred(output_size(A)) == (4,)
730+
# Test allocation-free interpolation with Vector{StaticArrays.SVector}
731+
u_s = [convert(SVector{length(u[1])}, i) for i in u]
732+
A_s = @inferred(QuadraticSpline(u_s, t; extrapolation = ExtrapolationType.Extension))
733+
for x in (-2.0, -0.5, 0.7, 2.0)
734+
@test A(x) == A_s(x)
735+
end
736+
@test A_s(0.7) isa SVector{length(u[1])}
737+
@test_nowarn test_allocs(A_s, 0.7)
683738

739+
# Test Vector{Matrix} interpolation
684740
u = [repeat(u[i], 1, 3) for i in 1:3]
685741
A = @inferred(QuadraticSpline(u, t; extrapolation = ExtrapolationType.Extension))
686742
@test A(-2.0) == P₁(-2.0) * ones(4, 3)
@@ -728,8 +784,9 @@ end
728784

729785
u_ = [0.0, 1.0, 3.0]' .* ones(4)
730786
u = [u_[:, i] for i in 1:size(u_, 2)]
787+
# Test Vector{Vector} interpolation
731788
@test @inferred(CubicSpline(u, t; extrapolation = ExtrapolationType.Extension)) isa
732-
CubicSpline broken=VERSION < v"1.11"
789+
CubicSpline
733790
A = CubicSpline(u, t; extrapolation = ExtrapolationType.Extension)
734791
for x in (-1.5, -0.5, -0.7)
735792
@test A(x) P₁(x) * ones(4)
@@ -739,9 +796,20 @@ end
739796
end
740797
@test @inferred(output_dim(A)) == 1
741798
@test @inferred(output_size(A)) == (4,)
799+
# Test allocation-free interpolation with StaticArrays
800+
u_s = [convert(SVector{length(first(u))}, i) for i in u]
801+
@test_broken @inferred(CubicSpline(
802+
u_s, t; extrapolation = ExtrapolationType.Extension)) isa CubicSpline
803+
A_s = CubicSpline(u_s, t; extrapolation = ExtrapolationType.Extension)
804+
for x in (-1.5, -0.5, -0.7)
805+
@test A(x) == A_s(x)
806+
end
807+
@test A_s(0) isa SVector{length(first(u))}
808+
@test_nowarn test_allocs(A_s, 0)
742809

810+
# Test Vector{Matrix} interpolation
743811
u = [repeat(u[i], 1, 3) for i in 1:3]
744-
@test_broken @inferred(CubicSpline(
812+
@test @inferred(CubicSpline(
745813
u, t; extrapolation = ExtrapolationType.Extension)) isa CubicSpline
746814
A = CubicSpline(u, t; extrapolation = ExtrapolationType.Extension)
747815
for x in (-1.5, -0.5, -0.7)
@@ -778,7 +846,7 @@ end
778846
0.0 cos(2t)]
779847
t = 0.1:0.1:1.0
780848
u3d = f3d.(t)
781-
@test_broken @inferred(CubicSpline(u3d, t)) isa CubicSpline
849+
@test @inferred(CubicSpline(u3d, t)) isa CubicSpline
782850
c = CubicSpline(u3d, t)
783851
t_test = 0.1:0.05:1.0
784852
u_test = reduce(hcat, c.(t_test))
@@ -952,8 +1020,18 @@ end
9521020
@testset "Vector of Vectors case" begin
9531021
u2 = [[u[i], u[i] + 1] for i in eachindex(u)]
9541022
du2 = [[du[i], du[i]] for i in eachindex(du)]
955-
A2 = CubicHermiteSpline(du2, u2, t)
1023+
A2 = CubicHermiteSpline(du2, u2, t; extrapolation = ExtrapolationType.Extension)
9561024
@test u2 A2.(t)
1025+
@test A2(100.0) repeat([10.106770], 2) + [0, 1] rtol=1e-5
1026+
@test A2(300.0) repeat([9.901542], 2) + [0, 1] rtol=1e-5
1027+
# Test allocation-free interpolation with Vector{StaticArrays.SVector}
1028+
u2_s = [convert(SVector{length(u2[1])}, i) for i in u2]
1029+
du2_s = [convert(SVector{length(du2[1])}, i) for i in du2]
1030+
A2_s = @inferred(CubicHermiteSpline(du2_s, u2_s, t; extrapolation = ExtrapolationType.Extension))
1031+
@test A2_s(100.0) == A2(100.0)
1032+
@test A2_s(300.0) == A2(300.0)
1033+
@test A2_s(0.7) isa SVector{length(u2[1])}
1034+
@test_nowarn test_allocs(A2_s, 0.7)
9571035
end
9581036
@testset "Vector of Matrices case" begin
9591037
u3 = [[u[i] u[i] + 1] for i in eachindex(u)]
@@ -1003,8 +1081,19 @@ end
10031081
u2 = [[u[i], u[i] + 1] for i in eachindex(u)]
10041082
du2 = [[du[i], du[i]] for i in eachindex(du)]
10051083
ddu2 = [[ddu[i], ddu[i]] for i in eachindex(ddu)]
1006-
A2 = QuinticHermiteSpline(ddu2, du2, u2, t)
1084+
A2 = QuinticHermiteSpline(ddu2, du2, u2, t; extrapolation = ExtrapolationType.Extension)
10071085
@test u2 A2.(t)
1086+
@test A2(100.0) repeat([10.107996], 2) + [0, 1] rtol=1e-5
1087+
@test A2(300.0) repeat([11.364162], 2) + [0, 1] rtol=1e-5
1088+
# Test allocation-free interpolation with Vector{StaticArrays.SVector}
1089+
u2_s = [convert(SVector{length(u2[1])}, i) for i in u2]
1090+
du2_s = [convert(SVector{length(du2[1])}, i) for i in du2]
1091+
ddu2_s = [convert(SVector{length(du2[1])}, i) for i in ddu2]
1092+
A2_s = @inferred(QuinticHermiteSpline(ddu2_s, du2_s, u2_s, t; extrapolation = ExtrapolationType.Extension))
1093+
@test A2_s(100.0) == A2(100.0)
1094+
@test A2_s(300.0) == A2(300.0)
1095+
@test A2_s(0.7) isa SVector{length(u2[1])}
1096+
@test_nowarn test_allocs(A2_s, 0.7)
10081097
end
10091098
@testset "Vector of Matrices case" begin
10101099
u3 = [[u[i] u[i] + 1] for i in eachindex(u)]

0 commit comments

Comments
 (0)