Skip to content
Merged
Show file tree
Hide file tree
Changes from 1 commit
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
1 change: 1 addition & 0 deletions Project.toml
Original file line number Diff line number Diff line change
Expand Up @@ -22,6 +22,7 @@ Reexport = "189a3867-3050-52da-a836-e630ba90ab69"
SparseArrays = "2f01184e-e22b-5df5-ae63-d93ebab69eaf"
SparseDiffTools = "47a9eef4-7e08-11e9-0b38-333d64bd3804"
StaticArrays = "90137ffa-7385-5640-81b9-e52037218182"
SuiteSparse = "4607b0f0-06f3-5cda-b6b1-a6196a1729e9"

[compat]
ArrayInterface = "1.1, 2.0"
Expand Down
2 changes: 1 addition & 1 deletion src/OrdinaryDiffEq.jl
Original file line number Diff line number Diff line change
Expand Up @@ -5,7 +5,7 @@ module OrdinaryDiffEq

using Logging

using MuladdMacro, SparseArrays
using MuladdMacro, SparseArrays, SuiteSparse

using LinearAlgebra

Expand Down
15 changes: 12 additions & 3 deletions src/derivative_utils.jl
Original file line number Diff line number Diff line change
Expand Up @@ -543,9 +543,18 @@ function build_J_W(alg,u,uprev,p,t,dt,f,uEltypeNoUnits,::Val{IIP}) where IIP
elseif u isa Number
u
else
LU{LinearAlgebra.lutype(uEltypeNoUnits)}(Matrix{uEltypeNoUnits}(undef, 0, 0),
Vector{LinearAlgebra.BlasInt}(undef, 0),
zero(LinearAlgebra.BlasInt))
# TODO: make this more general, maybe by running calc_W and use its returned value
if f.jac_prototype isa SparseMatrixCSC
SuiteSparse.UMFPACK.UmfpackLU(Ptr{Cvoid}(), Ptr{Cvoid}(), 1, 1,
f.jac_prototype.colptr[1:1],
f.jac_prototype.rowval[1:1],
f.jac_prototype.nzval[1:1],
0)
else
LU{LinearAlgebra.lutype(uEltypeNoUnits)}(Matrix{uEltypeNoUnits}(undef, 0, 0),
Vector{LinearAlgebra.BlasInt}(undef, 0),
zero(LinearAlgebra.BlasInt))
end
end
end # end W
end
Expand Down
3 changes: 1 addition & 2 deletions src/derivative_wrappers.jl
Original file line number Diff line number Diff line change
Expand Up @@ -41,7 +41,7 @@ function jacobian_autodiff(f, x::AbstractArray, odefun)
if DiffEqBase.has_colorvec(odefun)
colorvec = odefun.colorvec
sparsity = odefun.jac_prototype
jac_prototype = nothing
jac_prototype = odefun.jac_prototype
else
colorvec = 1:length(x)
sparsity = nothing
Expand Down Expand Up @@ -71,7 +71,6 @@ jacobian_finitediff(f, x, diff_type, dir, colorvec, sparsity, jac_prototype) =
jacobian_finitediff(f, x::AbstractArray, diff_type, dir, colorvec, sparsity, jac_prototype) =
(FiniteDiff.finite_difference_jacobian(f, x, diff_type, eltype(x), diff_type==Val{:forward} ? f(x) : similar(x),
dir = dir, colorvec = colorvec, sparsity = sparsity, jac_prototype = jac_prototype),_nfcount(maximum(colorvec),diff_type))

function jacobian(f, x, integrator)
alg = unwrap_alg(integrator, true)
local tmp
Expand Down