Skip to content

Commit 859c634

Browse files
committed
preallocated coordinate format vectors
1 parent 2f95d6b commit 859c634

File tree

3 files changed

+188
-23
lines changed

3 files changed

+188
-23
lines changed

src/julia_interface.jl

Lines changed: 75 additions & 23 deletions
Original file line numberDiff line numberDiff line change
@@ -2,6 +2,42 @@ export cons_coord, cons_coord!, consjac
22

33
using NLPModels, SparseArrays
44

5+
# Type conversion helpers for Issue #392: Preallocate more vectors
6+
"""
7+
prepare_input!(workspace::Vector{T}, x::AbstractVector{S}) where {T, S}
8+
9+
Prepare input vector `x` for use with CUTEst functions, using preallocated workspace
10+
to avoid allocations when type conversion is needed.
11+
12+
Returns the input vector directly if no conversion is needed, or the workspace vector
13+
with converted values if type conversion is required.
14+
"""
15+
@inline function prepare_input!(workspace::Vector{T}, x::AbstractVector{S}) where {T, S}
16+
if S === T && typeof(x) <: Vector{T}
17+
return x # No conversion needed
18+
else
19+
resize!(workspace, length(x))
20+
workspace .= x
21+
return workspace
22+
end
23+
end
24+
25+
"""
26+
prepare_output!(workspace::Vector{T}, target::AbstractVector{S}, source::AbstractVector{T}) where {T, S}
27+
28+
Prepare output by copying from source to target, using workspace for type conversion if needed.
29+
"""
30+
@inline function prepare_output!(workspace::Vector{T}, target::AbstractVector{S}, source::AbstractVector{T}) where {T, S}
31+
if S === T && typeof(target) <: Vector{T}
32+
target .= source
33+
else
34+
resize!(workspace, length(source))
35+
workspace .= source
36+
target .= workspace
37+
end
38+
return target
39+
end
40+
541
function NLPModels.obj(nlp::CUTEstModel{T}, x::StrideOneVector{T}) where {T}
642
@lencheck nlp.meta.nvar x
743
if nlp.meta.ncon > 0
@@ -40,10 +76,17 @@ end
4076

4177
function NLPModels.grad!(nlp::CUTEstModel{T}, x::AbstractVector, g::AbstractVector) where {T}
4278
@lencheck nlp.meta.nvar x g
43-
x_ = Vector{T}(x)
44-
g_ = Vector{T}(g)
45-
grad!(nlp, x_, g_)
46-
g .= g_
79+
80+
# Use type conversion helpers to avoid allocations (Issue #392)
81+
x_prepared = prepare_input!(nlp.input_workspace, x)
82+
83+
if typeof(g) <: Vector{T}
84+
grad!(nlp, x_prepared, g)
85+
else
86+
resize!(nlp.output_workspace, length(g))
87+
grad!(nlp, x_prepared, view(nlp.output_workspace, 1:length(g)))
88+
g .= view(nlp.output_workspace, 1:length(g))
89+
end
4790
end
4891

4992
function NLPModels.objcons!(
@@ -105,13 +148,15 @@ end
105148

106149
function NLPModels.objgrad!(nlp::CUTEstModel{T}, x::AbstractVector, g::StrideOneVector{T}) where {T}
107150
@lencheck nlp.meta.nvar x g
108-
objgrad!(nlp, Vector{T}(x), g)
151+
x_prepared = prepare_input!(nlp.input_workspace, x)
152+
objgrad!(nlp, x_prepared, g)
109153
end
110154

111155
function NLPModels.objgrad!(nlp::CUTEstModel{T}, x::AbstractVector, g::AbstractVector) where {T}
112156
@lencheck nlp.meta.nvar x g
157+
x_prepared = prepare_input!(nlp.input_workspace, x)
113158
gc = nlp.workspace_nvar
114-
f, _ = objgrad!(nlp, Vector{T}(x), gc)
159+
f, _ = objgrad!(nlp, x_prepared, gc)
115160
g .= gc
116161
return f, g
117162
end
@@ -179,15 +224,15 @@ function cons_coord!(
179224
@lencheck nlp.meta.nvar x
180225
@lencheck nlp.meta.ncon c
181226
@lencheck nlp.meta.nnzj rows cols vals
182-
rows_ = Vector{Cint}(undef, nlp.meta.nnzj)
183-
cols_ = Vector{Cint}(undef, nlp.meta.nnzj)
184-
vals_ = Vector{T}(undef, nlp.meta.nnzj)
185-
c_ = Vector{T}(undef, nlp.meta.ncon)
186-
cons_coord!(nlp, Vector{T}(x), c_, rows_, cols_, vals_)
187-
rows .= rows_
188-
cols .= cols_
189-
vals .= vals_
190-
c .= c_
227+
228+
# Use preallocated vectors instead of allocating new ones (Issue #392)
229+
cons_coord!(nlp, Vector{T}(x), nlp.cons_vals, nlp.jac_coord_rows, nlp.jac_coord_cols, nlp.jac_coord_vals)
230+
231+
# Copy results to output vectors
232+
rows .= nlp.jac_coord_rows
233+
cols .= nlp.jac_coord_cols
234+
vals .= nlp.jac_coord_vals
235+
c .= nlp.cons_vals
191236
return c, rows, cols, vals
192237
end
193238

@@ -209,11 +254,17 @@ Usage:
209254
"""
210255
function cons_coord(nlp::CUTEstModel{T}, x::StrideOneVector{T}) where {T}
211256
@lencheck nlp.meta.nvar x
212-
c = Vector{T}(undef, nlp.meta.ncon)
213-
rows = Vector{Cint}(undef, nlp.meta.nnzj)
214-
cols = Vector{Cint}(undef, nlp.meta.nnzj)
215-
vals = Vector{T}(undef, nlp.meta.nnzj)
216-
cons_coord!(nlp, x, c, rows, cols, vals)
257+
258+
# Use preallocated vectors to avoid allocations (Issue #392)
259+
cons_coord!(nlp, x, nlp.cons_vals, nlp.jac_coord_rows, nlp.jac_coord_cols, nlp.jac_coord_vals)
260+
261+
# Return copies of the results to maintain API compatibility
262+
c = copy(nlp.cons_vals)
263+
rows = copy(nlp.jac_coord_rows)
264+
cols = copy(nlp.jac_coord_cols)
265+
vals = copy(nlp.jac_coord_vals)
266+
267+
return c, rows, cols, vals
217268
end
218269

219270
function cons_coord(nlp::CUTEstModel{T}, x::AbstractVector) where {T}
@@ -630,9 +681,10 @@ function NLPModels.hess_coord!(
630681
@lencheck nlp.meta.nvar x
631682
@lencheck nlp.meta.ncon y
632683
@lencheck nlp.meta.nnzh vals
633-
vals_ = Vector{T}(undef, nlp.meta.nnzh)
634-
NLPModels.hess_coord!(nlp, Vector{T}(x), convert(Vector{T}, y), vals_, obj_weight = obj_weight)
635-
vals .= vals_
684+
685+
# Use preallocated vector instead of allocating (Issue #392)
686+
NLPModels.hess_coord!(nlp, Vector{T}(x), convert(Vector{T}, y), nlp.hess_coord_vals, obj_weight = obj_weight)
687+
vals .= nlp.hess_coord_vals
636688
return vals
637689
end
638690

src/model.jl

Lines changed: 37 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -15,6 +15,20 @@ mutable struct CUTEstModel{T} <: AbstractNLPModel{T, Vector{T}}
1515
workspace_nvar::Vector{T}
1616
workspace_ncon::Vector{T}
1717

18+
# Preallocated coordinate format vectors
19+
jac_coord_rows::Vector{Cint} # nnzj elements for Jacobian row indices
20+
jac_coord_cols::Vector{Cint} # nnzj elements for Jacobian column indices
21+
jac_coord_vals::Vector{T} # nnzj elements for Jacobian values
22+
hess_coord_vals::Vector{T} # nnzh elements for Hessian values
23+
24+
# Preallocated constraint evaluation vectors
25+
cons_vals::Vector{T} # ncon elements for constraint values
26+
cons_nln_vals::Vector{T} # nnln elements for nonlinear constraints subset
27+
28+
# Type conversion workspace vectors
29+
input_workspace::Vector{T} # nvar elements for input conversion
30+
output_workspace::Vector{T} # max(nvar, ncon) elements for output conversion
31+
1832
Jval::Vector{T}
1933
Jvar::Vector{Cint}
2034

@@ -297,6 +311,21 @@ function CUTEstModel{T}(
297311
workspace_nvar = Vector{T}(undef, nvar)
298312
workspace_ncon = Vector{T}(undef, ncon)
299313

314+
# Preallocate new coordinate format vectors (Issue #392)
315+
jac_coord_rows = Vector{Cint}(undef, nnzj)
316+
jac_coord_cols = Vector{Cint}(undef, nnzj)
317+
jac_coord_vals = Vector{T}(undef, nnzj)
318+
hess_coord_vals = Vector{T}(undef, nnzh)
319+
320+
# Preallocate constraint evaluation vectors
321+
cons_vals = Vector{T}(undef, ncon)
322+
nnln = count(.!linear) # Number of nonlinear constraints
323+
cons_nln_vals = Vector{T}(undef, nnln)
324+
325+
# Preallocate type conversion workspace vectors
326+
input_workspace = Vector{T}(undef, nvar)
327+
output_workspace = Vector{T}(undef, max(nvar, ncon))
328+
300329
fclose(T, libsif, funit, status)
301330
cutest_error(status[])
302331

@@ -330,6 +359,14 @@ function CUTEstModel{T}(
330359
clinvals,
331360
workspace_nvar,
332361
workspace_ncon,
362+
jac_coord_rows,
363+
jac_coord_cols,
364+
jac_coord_vals,
365+
hess_coord_vals,
366+
cons_vals,
367+
cons_nln_vals,
368+
input_workspace,
369+
output_workspace,
333370
Jval,
334371
Jvar,
335372
libsif,

test/test_comprehensive.jl

Lines changed: 76 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,76 @@
1+
using CUTEst
2+
using BenchmarkTools
3+
using LinearAlgebra
4+
using NLPModels
5+
6+
# Test with an unconstrained problem first
7+
println("Testing with unconstrained problem (ROSENBR):")
8+
nlp1 = CUTEstModel{Float64}("ROSENBR")
9+
x0 = nlp1.meta.x0
10+
11+
# Test hess_coord allocations
12+
println("✓ Problem loaded: $(nlp1.meta.name) ($(nlp1.meta.nvar) vars, $(nlp1.meta.ncon) cons)")
13+
vals_h = Vector{Float64}(undef, nlp1.meta.nnzh)
14+
hess_coord!(nlp1, x0, vals_h)
15+
16+
b_hess = @benchmark hess_coord!($nlp1, $x0, $vals_h) samples=100 evals=5
17+
println("✓ hess_coord!: $(b_hess.allocs) allocations, $(BenchmarkTools.prettytime(median(b_hess).time))")
18+
19+
finalize(nlp1)
20+
21+
# Test with a constrained problem
22+
println("\nTesting with constrained problem (HS6):")
23+
nlp2 = CUTEstModel{Float64}("HS6")
24+
x0 = nlp2.meta.x0
25+
26+
println("✓ Problem loaded: $(nlp2.meta.name) ($(nlp2.meta.nvar) vars, $(nlp2.meta.ncon) cons)")
27+
println(" nnzj: $(nlp2.meta.nnzj), nnzh: $(nlp2.meta.nnzh)")
28+
29+
# Test basic functionality
30+
f = obj(nlp2, x0)
31+
g = similar(x0)
32+
grad!(nlp2, x0, g)
33+
println("✓ obj/grad: f = $(round(f, digits=6)), ||g|| = $(round(norm(g), digits=6))")
34+
35+
# Test constraint functions
36+
c = Vector{Float64}(undef, nlp2.meta.ncon)
37+
cons!(nlp2, x0, c)
38+
println("✓ constraints: ||c|| = $(round(norm(c), digits=6))")
39+
40+
# Test cons_coord - this should show major allocation improvements
41+
println("\nTesting cons_coord allocation improvements:")
42+
c1, rows1, cols1, vals1 = cons_coord(nlp2, x0)
43+
44+
# Benchmark cons_coord
45+
b_cons = @benchmark cons_coord($nlp2, $x0) samples=100 evals=5
46+
println("✓ cons_coord: $(b_cons.allocs) allocations, $(BenchmarkTools.prettytime(median(b_cons).time))")
47+
println(" Memory: $(BenchmarkTools.prettymemory(median(b_cons).memory))")
48+
println(" Returned $(length(vals1)) Jacobian elements")
49+
50+
# Test cons_coord!
51+
rows = Vector{Cint}(undef, nlp2.meta.nnzj)
52+
cols = Vector{Cint}(undef, nlp2.meta.nnzj)
53+
vals = Vector{Float64}(undef, nlp2.meta.nnzj)
54+
c_out = Vector{Float64}(undef, nlp2.meta.ncon)
55+
56+
b_cons_inplace = @benchmark cons_coord!($nlp2, $x0, $c_out, $rows, $cols, $vals) samples=100 evals=5
57+
println("✓ cons_coord!: $(b_cons_inplace.allocs) allocations, $(BenchmarkTools.prettytime(median(b_cons_inplace).time))")
58+
59+
# Test type conversion
60+
println("\nTesting type conversion improvements:")
61+
x0_f32 = Float32.(x0)
62+
g_f32 = Vector{Float32}(undef, nlp2.meta.nvar)
63+
64+
b_grad_conv = @benchmark grad!($nlp2, $x0_f32, $g_f32) samples=100 evals=5
65+
println("✓ grad! with Float32->Float64 conversion: $(b_grad_conv.allocs) allocations")
66+
println(" Time: $(BenchmarkTools.prettytime(median(b_grad_conv).time))")
67+
68+
# Test hess_coord with constraints
69+
vals_h2 = Vector{Float64}(undef, nlp2.meta.nnzh)
70+
y = zeros(nlp2.meta.ncon)
71+
hess_coord!(nlp2, x0, y, vals_h2)
72+
73+
b_hess2 = @benchmark hess_coord!($nlp2, $x0, $y, $vals_h2) samples=100 evals=5
74+
println("✓ hess_coord! (constrained): $(b_hess2.allocs) allocations, $(BenchmarkTools.prettytime(median(b_hess2).time))")
75+
76+
finalize(nlp2)

0 commit comments

Comments
 (0)