Skip to content
This repository was archived by the owner on Aug 22, 2025. It is now read-only.

Commit b1a3244

Browse files
author
Stuart Daines
committed
Add fast path for forwarddiff_color_jacobian! with sparse J
Partial fix for #138 Bodges in a fast path for the case where J and sparsity are are both SparseMatrixCSC and have the same number of columns and stored values.
1 parent e798411 commit b1a3244

1 file changed

Lines changed: 32 additions & 1 deletion

File tree

src/differentiation/compute_jacobian_ad.jl

Lines changed: 32 additions & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -264,8 +264,20 @@ function forwarddiff_color_jacobian!(J::AbstractMatrix{<:Number},
264264
chunksize = jac_cache.chunksize
265265
color_i = 1
266266
maxcolor = maximum(colorvec)
267+
267268
fill!(J, zero(eltype(J)))
268269

270+
# fast path if J and sparsity are both SparseMatrixCSC and have the same number of columns and stored values
271+
sparseCSC_common = (J isa SparseMatrixCSC &&
272+
sparsity isa SparseMatrixCSC &&
273+
length(J.colptr) == length(sparsity.colptr) &&
274+
length(J.nzval) == length(sparsity.nzval))
275+
276+
if sparseCSC_common
277+
J.colptr .= sparsity.colptr
278+
J.rowval .= sparsity.rowval
279+
end
280+
269281
if FiniteDiff._use_findstructralnz(sparsity)
270282
rows_index, cols_index = ArrayInterface.findstructralnz(sparsity)
271283
else
@@ -287,9 +299,14 @@ function forwarddiff_color_jacobian!(J::AbstractMatrix{<:Number},
287299
if !(sparsity isa Nothing)
288300
for j in 1:chunksize
289301
dx .= partials.(fx, j)
302+
290303
if ArrayInterface.fast_scalar_indexing(dx)
291304
#dx is implicitly used in vecdx
292-
FiniteDiff._colorediteration!(J,sparsity,rows_index,cols_index,vecdx,colorvec,color_i,ncols)
305+
if sparseCSC_common
306+
_colorediteration!(J,vecdx,colorvec,color_i,ncols)
307+
else
308+
FiniteDiff._colorediteration!(J,sparsity,rows_index,cols_index,vecdx,colorvec,color_i,ncols)
309+
end
293310
else
294311
#=
295312
J.nzval[rows_index] .+= (colorvec[cols_index] .== color_i) .* dx[rows_index]
@@ -316,3 +333,17 @@ function forwarddiff_color_jacobian!(J::AbstractMatrix{<:Number},
316333
end
317334
return J
318335
end
336+
337+
# fast version of FiniteDiff._colorediteration! for the case where J and sparsity have the same sparsity pattern
338+
@inline function _colorediteration!(Jsparsity::SparseMatrixCSC,vfx,colorvec,color_i,ncols)
339+
@inbounds for col_index in 1:ncols
340+
if colorvec[col_index] == color_i
341+
@inbounds for spidx in nzrange(Jsparsity, col_index)
342+
row_index = Jsparsity.rowval[spidx]
343+
Jsparsity.nzval[spidx]=vfx[row_index]
344+
end
345+
end
346+
end
347+
end
348+
349+

0 commit comments

Comments
 (0)