Skip to content

Commit 30261fc

Browse files
committed
Don't create a reshaped array when compute a Jacobian of a StaticVector
valued function of size one.
1 parent 310da8c commit 30261fc

File tree

4 files changed

+45
-30
lines changed

4 files changed

+45
-30
lines changed

Project.toml

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -19,8 +19,8 @@ julia = "1.2"
1919
BandedMatrices = "aae01518-5342-5314-be14-df237901396f"
2020
BlockBandedMatrices = "ffab5731-97b5-5995-9138-79e8c1846df0"
2121
Pkg = "44cfe95a-1eb2-52ea-b672-e2afdf69b78f"
22-
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
2322
SafeTestsets = "1bc83da4-3b8d-516f-aca4-4fe02f6d838f"
23+
Test = "8dfed614-e22c-5e08-85e1-65c5234f0b40"
2424

2525
[targets]
2626
test = ["Test", "BlockBandedMatrices", "BandedMatrices", "Pkg", "SafeTestsets"]

src/FiniteDiff.jl

Lines changed: 4 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -8,6 +8,10 @@ import Base: resize!
88
_vec(x) = vec(x)
99
_vec(x::Number) = x
1010

11+
_mat(x::AbstractMatrix) = x
12+
_mat(x::StaticVector) = reshape(x, (axes(x, 1), SOneTo(1)))
13+
_mat(x::AbstractVector) = reshape(x, (axes(x, 1), OneTo(1)))
14+
1115
include("iteration_utils.jl")
1216
include("epsilons.jl")
1317
include("derivatives.jl")

src/jacobians.jl

Lines changed: 15 additions & 21 deletions
Original file line numberDiff line numberDiff line change
@@ -121,7 +121,7 @@ end
121121

122122
function _make_Ji(::AbstractArray, xtype, dx, color_i, nrows, ncols)
123123
Ji = mapreduce(i -> i==color_i ? dx : zero(dx), hcat, 1:ncols)
124-
size(Ji)!=(nrows, ncols) ? reshape(Ji,(nrows,ncols)) : Ji #branch when size(dx) == (1,) => size(Ji) == (1,) while size(J) == (1,1)
124+
size(Ji) != (nrows, ncols) ? reshape(Ji, (nrows, ncols)) : Ji #branch when size(dx) == (1,) => size(Ji) == (1,) while size(J) == (1,1)
125125
end
126126

127127
function finite_difference_jacobian(f, x,
@@ -186,17 +186,15 @@ function finite_difference_jacobian(
186186
x_save = ArrayInterface.allowed_getindex(vecx, i)
187187
epsilon = compute_epsilon(Val{:forward}, x_save, relstep, absstep, dir)
188188
_vecx1 = Base.setindex(vecx, x_save+epsilon, i)
189-
_x1 = reshape(_vecx1,size(x))
189+
_x1 = reshape(_vecx1, axes(x))
190190
vecfx1 = _vec(f(_x1))
191191
dx = (vecfx1-vecfx) / epsilon
192192
return dx
193193
end
194194

195195
if jac_prototype isa Nothing && sparsity isa Nothing
196196
J = mapreduce(calculate_Ji_forward, hcat, 1:maximum(colorvec))
197-
if maximum(colorvec) == 1
198-
J = reshape(J, 1, 1)
199-
end
197+
J = _mat(J)
200198
else
201199
@inbounds for color_i 1:maximum(colorvec)
202200
if sparsity isa Nothing
@@ -206,7 +204,7 @@ function finite_difference_jacobian(
206204
tmp = norm(vecx .* (colorvec .== color_i))
207205
epsilon = compute_epsilon(Val{:forward}, sqrt(tmp), relstep, absstep, dir)
208206
_vecx = @. vecx + epsilon * (colorvec == color_i)
209-
_x = reshape(_vecx,size(x))
207+
_x = reshape(_vecx, axes(x))
210208
vecfx1 = _vec(f(_x))
211209
dx = (vecfx1-vecfx)/epsilon
212210
Ji = _make_Ji(J,rows_index,cols_index,dx,colorvec,color_i,nrows,ncols)
@@ -222,8 +220,8 @@ function finite_difference_jacobian(
222220
epsilon = compute_epsilon(Val{:forward}, x1_save, relstep, absstep, dir)
223221
_vecx1 = Base.setindex(vecx1,x1_save+epsilon,i)
224222
_vecx = Base.setindex(vecx,x_save-epsilon,i)
225-
_x1 = reshape(_vecx1,size(x))
226-
_x = reshape(_vecx,size(x))
223+
_x1 = reshape(_vecx1, axes(x))
224+
_x = reshape(_vecx, axes(x))
227225
vecfx1 = _vec(f(_x1))
228226
vecfx = _vec(f(_x))
229227
dx = (vecfx1-vecfx)/(2epsilon)
@@ -232,9 +230,7 @@ function finite_difference_jacobian(
232230

233231
if jac_prototype isa Nothing && sparsity isa Nothing
234232
J = mapreduce(calculate_Ji_central, hcat, 1:maximum(colorvec))
235-
if maximum(colorvec) == 1
236-
J = reshape(J, 1, 1)
237-
end
233+
J = _mat(J)
238234
else
239235
@inbounds for color_i 1:maximum(colorvec)
240236
if sparsity isa Nothing
@@ -245,8 +241,8 @@ function finite_difference_jacobian(
245241
epsilon = compute_epsilon(Val{:forward}, sqrt(tmp), relstep, absstep, dir)
246242
_vecx1 = @. vecx1 + epsilon * (colorvec == color_i)
247243
_vecx = @. vecx - epsilon * (colorvec == color_i)
248-
_x1 = reshape(_vecx1,size(x))
249-
_x = reshape(_vecx, size(x))
244+
_x1 = reshape(_vecx1, axes(x))
245+
_x = reshape(_vecx, axes(x))
250246
vecfx1 = _vec(f(_x1))
251247
vecfx = _vec(f(_x))
252248
dx = (vecfx1-vecfx)/(2epsilon)
@@ -257,29 +253,27 @@ function finite_difference_jacobian(
257253
end
258254
elseif fdtype == Val{:complex} && returntype <: Real
259255
epsilon = eps(eltype(x))
260-
256+
261257
function calculate_Ji_complex(i)
262258
x_save = ArrayInterface.allowed_getindex(vecx,i)
263259
_vecx = Base.setindex(complex.(vecx),x_save+im*epsilon,i)
264-
_x = reshape(_vecx,size(x))
260+
_x = reshape(_vecx, axes(x))
265261
vecfx = _vec(f(_x))
266262
dx = imag(vecfx)/epsilon
267263
return dx
268264
end
269-
265+
270266
if jac_prototype isa Nothing && sparsity isa Nothing
271267
J = mapreduce(calculate_Ji_complex, hcat, 1:maximum(colorvec))
272-
if maximum(colorvec) == 1
273-
J = reshape(J, 1, 1)
274-
end
268+
J = _mat(J)
275269
else
276270
@inbounds for color_i 1:maximum(colorvec)
277271
if sparsity isa Nothing
278272
dx = calculate_Ji_complex(color_i)
279273
J = J + _make_Ji(J, eltype(x), dx, color_i, nrows, ncols)
280274
else
281275
_vecx = @. vecx + im * epsilon * (colorvec == color_i)
282-
_x = reshape(_vecx,size(x))
276+
_x = reshape(_vecx, axes(x))
283277
vecfx = _vec(f(_x))
284278
dx = imag(vecfx)/epsilon
285279
Ji = _make_Ji(J,rows_index,cols_index,dx,colorvec,color_i,nrows,ncols)
@@ -332,7 +326,7 @@ function finite_difference_jacobian!(
332326
dir = true) where {T1,T2,T3,cType,sType,fdtype,returntype}
333327

334328
m, n = size(J)
335-
_color = reshape(colorvec,size(x)...)
329+
_color = reshape(colorvec, axes(x)...)
336330

337331
x1, fx, fx1 = cache.x1, cache.fx, cache.fx1
338332
copyto!(x1, x)

test/out_of_place_tests.jl

Lines changed: 25 additions & 8 deletions
Original file line numberDiff line numberDiff line change
@@ -16,28 +16,45 @@ function second_derivative_stencil(N)
1616
end
1717

1818
x = @SVector ones(30)
19-
J = FiniteDiff.finite_difference_jacobian(f,x, Val{:forward}, eltype(x))
19+
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:forward}, eltype(x))
2020
@test J second_derivative_stencil(30)
2121

22-
J = FiniteDiff.finite_difference_jacobian(f,x, Val{:central}, eltype(x))
22+
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:central}, eltype(x))
2323
@test J second_derivative_stencil(30)
2424

25-
J = FiniteDiff.finite_difference_jacobian(f,x, Val{:complex}, eltype(x))
25+
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:complex}, eltype(x))
2626
@test J second_derivative_stencil(30)
2727

2828
spJ = sparse(second_derivative_stencil(30))
29-
J = FiniteDiff.finite_difference_jacobian(f,x, Val{:forward}, eltype(x),jac_prototype=spJ)
29+
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:forward}, eltype(x), jac_prototype=spJ)
3030
@test J second_derivative_stencil(30)
3131
@test typeof(J) == typeof(spJ)
32-
J = FiniteDiff.finite_difference_jacobian(f,x, Val{:forward}, eltype(x),colorvec=repeat(1:3,10),sparsity=spJ,jac_prototype=spJ)
32+
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:forward}, eltype(x),
33+
colorvec=SVector{30}(repeat(1:3, 10)), sparsity=spJ, jac_prototype=spJ)
3334
@test J second_derivative_stencil(30)
3435
@test typeof(J) == typeof(spJ)
36+
3537
#1x1 SVector test
3638
x = SVector{1}([1.])
3739
f(x) = x
3840
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:forward}, eltype(x))
39-
@test J SMatrix{1,1}([1.])
41+
@test J[1, 1] 1.0
42+
@test J isa SMatrix{1,1}
43+
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:central}, eltype(x))
44+
@test J[1, 1] 1.0
45+
@test J isa SMatrix{1,1}
46+
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:complex}, eltype(x))
47+
@test J[1, 1] 1.0
48+
@test J isa SMatrix{1,1}
49+
50+
x = SVector{1}([1.])
51+
f(x) = vcat(x, x)
52+
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:forward}, eltype(x))
53+
@test J fill(1.0, 2, 1)
54+
@test J isa SMatrix{2,1}
4055
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:central}, eltype(x))
41-
@test J SMatrix{1,1}([1.])
56+
@test J fill(1.0, 2, 1)
57+
@test J isa SMatrix{2,1}
4258
J = FiniteDiff.finite_difference_jacobian(f, x, Val{:complex}, eltype(x))
43-
@test J SMatrix{1,1}([1.])
59+
@test J fill(1.0, 2, 1)
60+
@test J isa SMatrix{2,1}

0 commit comments

Comments
 (0)