Differentiation of Simple ODE Benchmarks
From the paper A Comparison of Automatic Differentiation and Continuous Sensitivity Analysis for Derivatives of Differential Equation Solutions
using ParameterizedFunctions, OrdinaryDiffEq, LinearAlgebra, StaticArrays
using SciMLSensitivity, ForwardDiff, FiniteDiff, ReverseDiff, BenchmarkTools, Test
using DataFrames, PrettyTables, Markdown
tols = (abstol=1e-5, reltol=1e-7)
(abstol = 1.0e-5, reltol = 1.0e-7)
Define the Test ODEs
function lvdf(du, u, p, t)
a,b,c = p
x, y = u
du[1] = a*x - b*x*y
du[2] = -c*y + x*y
nothing
end
function lvcom_df(du, u, p, t)
a,b,c = p
x, y, s1, s2, s3, s4, s5, s6 = u
du[1] = a*x - b*x*y
du[2] = -c*y + x*y
#####################
# [a-by -bx]
# J = [ ]
# [y x-c]
#####################
J = @SMatrix [a-b*y -b*x
y x-c]
JS = J*@SMatrix[s1 s3 s5
s2 s4 s6]
G = @SMatrix [x -x*y 0
0 0 -y]
du[3:end] .= vec(JS+G)
nothing
end
lvdf_with_jacobian = ODEFunction{true, SciMLBase.FullSpecialize}(lvdf, jac=(J,u,p,t)->begin
a,b,c = p
x, y = u
J[1] = a-b*y
J[2] = y
J[3] = -b*x
J[4] = x-c
nothing
end)
u0 = [1.,1.]; tspan = (0., 10.); p = [1.5,1.0,3.0]; lvcom_u0 = [u0...;zeros(6)]
lvprob = ODEProblem{true, SciMLBase.FullSpecialize}(lvcom_df, lvcom_u0, tspan, p)
ODEProblem with uType Vector{Float64} and tType Float64. In-place: true
timespan: (0.0, 10.0)
u0: 8-element Vector{Float64}:
1.0
1.0
0.0
0.0
0.0
0.0
0.0
0.0
pkpdf = @ode_def begin
dEv = -Ka1*Ev
dCent = Ka1*Ev - (CL+Vmax/(Km+(Cent/Vc))+Q)*(Cent/Vc) + Q*(Periph/Vp) - Q2*(Cent/Vc) + Q2*(Periph2/Vp2)
dPeriph = Q*(Cent/Vc) - Q*(Periph/Vp)
dPeriph2 = Q2*(Cent/Vc) - Q2*(Periph2/Vp2)
dResp = Kin*(1-(IMAX*(Cent/Vc)^γ/(IC50^γ+(Cent/Vc)^γ))) - Kout*Resp
end Ka1 CL Vc Q Vp Kin Kout IC50 IMAX γ Vmax Km Q2 Vp2
pkpdp = [
1, # Ka1 Absorption rate constant 1 (1/time)
1, # CL Clearance (volume/time)
20, # Vc Central volume (volume)
2, # Q Inter-compartmental clearance (volume/time)
10, # Vp Peripheral volume of distribution (volume)
10, # Kin Response in rate constant (1/time)
2, # Kout Response out rate constant (1/time)
2, # IC50 Concentration for 50% of max inhibition (mass/volume)
1, # IMAX Maximum inhibition
1, # γ Emax model sigmoidicity
0, # Vmax Maximum reaction velocity (mass/time)
2, # Km Michaelis constant (mass/volume)
0.5, # Q2 Inter-compartmental clearance2 (volume/time)
100 # Vp2 Peripheral2 volume of distribution (volume)
];
pkpdu0 = [100, eps(), eps(), eps(), 5.] # exact zero in the initial condition triggers NaN in Jacobian
#pkpdu0 = ones(5)
pkpdcondition = function (u,t,integrator)
t in 0:24:240
end
pkpdaffect! = function (integrator)
integrator.u[1] += 100
end
pkpdcb = DiscreteCallback(pkpdcondition, pkpdaffect!, save_positions=(false, true))
pkpdtspan = (0.,240.)
pkpdprob = ODEProblem{true, SciMLBase.FullSpecialize}(pkpdf.f, pkpdu0, pkpdtspan, pkpdp)
pkpdfcomp = let pkpdf=pkpdf, J=zeros(5,5), JP=zeros(5,14), tmpdu=zeros(5,14)
function (du, u, p, t)
pkpdf.f(@view(du[:, 1]), u, p, t)
pkpdf.jac(J, u, p, t)
pkpdf.paramjac(JP, u, p, t)
mul!(tmpdu, J, @view(u[:, 2:end]))
du[:, 2:end] .= tmpdu .+ JP
nothing
end
end
pkpdcompprob = ODEProblem{true, SciMLBase.FullSpecialize}(pkpdfcomp, hcat(pkpdprob.u0,zeros(5,14)), pkpdprob.tspan, pkpdprob.p)
ODEProblem with uType Matrix{Float64} and tType Float64. In-place: true
timespan: (0.0, 240.0)
u0: 5×15 Matrix{Float64}:
100.0 0.0 0.0 0.0 0.0 0.0 … 0.0 0.0 0.0 0.0 0.0 0.0 0
.0
2.22045e-16 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0
.0
2.22045e-16 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0
.0
2.22045e-16 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0
.0
5.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0 0
.0
pollution = @ode_def begin
dy1 = -k1 *y1-k10*y11*y1-k14*y1*y6-k23*y1*y4-k24*y19*y1+
k2 *y2*y4+k3 *y5*y2+k9 *y11*y2+k11*y13+k12*y10*y2+k22*y19+k25*y20
dy2 = -k2 *y2*y4-k3 *y5*y2-k9 *y11*y2-k12*y10*y2+k1 *y1+k21*y19
dy3 = -k15*y3+k1 *y1+k17*y4+k19*y16+k22*y19
dy4 = -k2 *y2*y4-k16*y4-k17*y4-k23*y1*y4+k15*y3
dy5 = -k3 *y5*y2+k4 *y7+k4 *y7+k6 *y7*y6+k7 *y9+k13*y14+k20*y17*y6
dy6 = -k6 *y7*y6-k8 *y9*y6-k14*y1*y6-k20*y17*y6+k3 *y5*y2+k18*y16+k18*y16
dy7 = -k4 *y7-k5 *y7-k6 *y7*y6+k13*y14
dy8 = k4 *y7+k5 *y7+k6 *y7*y6+k7 *y9
dy9 = -k7 *y9-k8 *y9*y6
dy10 = -k12*y10*y2+k7 *y9+k9 *y11*y2
dy11 = -k9 *y11*y2-k10*y11*y1+k8 *y9*y6+k11*y13
dy12 = k9 *y11*y2
dy13 = -k11*y13+k10*y11*y1
dy14 = -k13*y14+k12*y10*y2
dy15 = k14*y1*y6
dy16 = -k18*y16-k19*y16+k16*y4
dy17 = -k20*y17*y6
dy18 = k20*y17*y6
dy19 = -k21*y19-k22*y19-k24*y19*y1+k23*y1*y4+k25*y20
dy20 = -k25*y20+k24*y19*y1
end k1 k2 k3 k4 k5 k6 k7 k8 k9 k10 k11 k12 k13 k14 k15 k16 k17 k18 k19 k20 k21 k22 k23 k24 k25
function make_pollution()
comp = let pollution = pollution, J = zeros(20, 20), JP = zeros(20, 25), tmpdu = zeros(20,25), tmpu = zeros(20,25)
function comp(du, u, p, t)
tmpu .= @view(u[:, 2:26])
pollution(@view(du[:, 1]), u, p, t)
pollution.jac(J,u,p,t)
pollution.paramjac(JP,u,p,t)
mul!(tmpdu, J, tmpu)
du[:,2:26] .= tmpdu .+ JP
nothing
end
end
u0 = zeros(20)
p = [.35e0, .266e2, .123e5, .86e-3, .82e-3, .15e5, .13e-3, .24e5, .165e5, .9e4, .22e-1, .12e5, .188e1, .163e5, .48e7, .35e-3, .175e-1, .1e9, .444e12, .124e4, .21e1, .578e1, .474e-1, .178e4, .312e1]
u0[2] = 0.2
u0[4] = 0.04
u0[7] = 0.1
u0[8] = 0.3
u0[9] = 0.01
u0[17] = 0.007
compu0 = zeros(20, 26)
compu0[1:20] .= u0
comp, u0, p, compu0
end
make_pollution (generic function with 1 method)
function makebrusselator(N=8)
xyd_brusselator = range(0,stop=1,length=N)
function limit(a, N)
if a == N+1
return 1
elseif a == 0
return N
else
return a
end
end
brusselator_f(x, y, t) = ifelse((((x-0.3)^2 + (y-0.6)^2) <= 0.1^2) &&
(t >= 1.1), 5., 0.)
brusselator_2d_loop = let N=N, xyd=xyd_brusselator, dx=step(xyd_brusselator)
function brusselator_2d_loop(du, u, p, t)
@inbounds begin
ii1 = N^2
ii2 = ii1+N^2
ii3 = ii2+2(N^2)
A = @view p[1:ii1]
B = @view p[ii1+1:ii2]
α = @view p[ii2+1:ii3]
II = LinearIndices((N, N, 2))
for I in CartesianIndices((N, N))
x = xyd[I[1]]
y = xyd[I[2]]
i = I[1]
j = I[2]
ip1 = limit(i+1, N); im1 = limit(i-1, N)
jp1 = limit(j+1, N); jm1 = limit(j-1, N)
du[II[i,j,1]] = α[II[i,j,1]]*(u[II[im1,j,1]] + u[II[ip1,j,1]] + u[II[i,jp1,1]] + u[II[i,jm1,1]] - 4u[II[i,j,1]])/dx^2 +
B[II[i,j,1]] + u[II[i,j,1]]^2*u[II[i,j,2]] - (A[II[i,j,1]] + 1)*u[II[i,j,1]] + brusselator_f(x, y, t)
end
for I in CartesianIndices((N, N))
i = I[1]
j = I[2]
ip1 = limit(i+1, N)
im1 = limit(i-1, N)
jp1 = limit(j+1, N)
jm1 = limit(j-1, N)
du[II[i,j,2]] = α[II[i,j,2]]*(u[II[im1,j,2]] + u[II[ip1,j,2]] + u[II[i,jp1,2]] + u[II[i,jm1,2]] - 4u[II[i,j,2]])/dx^2 +
A[II[i,j,1]]*u[II[i,j,1]] - u[II[i,j,1]]^2*u[II[i,j,2]]
end
return nothing
end
end
end
function init_brusselator_2d(xyd)
N = length(xyd)
u = zeros(N, N, 2)
for I in CartesianIndices((N, N))
x = xyd[I[1]]
y = xyd[I[2]]
u[I,1] = 22*(y*(1-y))^(3/2)
u[I,2] = 27*(x*(1-x))^(3/2)
end
vec(u)
end
dx = step(xyd_brusselator)
e1 = ones(N-1)
off = N-1
e4 = ones(N-off)
T = diagm(0=>-2ones(N), -1=>e1, 1=>e1, off=>e4, -off=>e4) ./ dx^2
Ie = Matrix{Float64}(I, N, N)
# A + df/du
Op = kron(Ie, T) + kron(T, Ie)
brusselator_jac = let N=N
(J,a,p,t) -> begin
ii1 = N^2
ii2 = ii1+N^2
ii3 = ii2+2(N^2)
A = @view p[1:ii1]
B = @view p[ii1+1:ii2]
α = @view p[ii2+1:ii3]
u = @view a[1:end÷2]
v = @view a[end÷2+1:end]
N2 = length(a)÷2
α1 = @view α[1:end÷2]
α2 = @view α[end÷2+1:end]
fill!(J, 0)
J[1:N2, 1:N2] .= α1.*Op
J[N2+1:end, N2+1:end] .= α2.*Op
J1 = @view J[1:N2, 1:N2]
J2 = @view J[N2+1:end, 1:N2]
J3 = @view J[1:N2, N2+1:end]
J4 = @view J[N2+1:end, N2+1:end]
J1[diagind(J1)] .+= @. 2u*v-(A+1)
J2[diagind(J2)] .= @. A-2u*v
J3[diagind(J3)] .= @. u^2
J4[diagind(J4)] .+= @. -u^2
nothing
end
end
Jmat = zeros(2N*N, 2N*N)
dp = zeros(2N*N, 4N*N)
brusselator_comp = let N=N, xyd=xyd_brusselator, dx=step(xyd_brusselator), Jmat=Jmat, dp=dp, brusselator_jac=brusselator_jac
function brusselator_comp(dus, us, p, t)
@inbounds begin
ii1 = N^2
ii2 = ii1+N^2
ii3 = ii2+2(N^2)
@views u, s = us[1:ii2], us[ii2+1:end]
du = @view dus[1:ii2]
ds = @view dus[ii2+1:end]
fill!(dp, 0)
A = @view p[1:ii1]
B = @view p[ii1+1:ii2]
α = @view p[ii2+1:ii3]
dfdα = @view dp[:, ii2+1:ii3]
diagind(dfdα)
for i in 1:ii1
dp[i, ii1+i] = 1
end
II = LinearIndices((N, N, 2))
uu = @view u[1:end÷2]
for i in eachindex(uu)
dp[i, i] = -uu[i]
dp[i+ii1, i] = uu[i]
end
for I in CartesianIndices((N, N))
x = xyd[I[1]]
y = xyd[I[2]]
i = I[1]
j = I[2]
ip1 = limit(i+1, N); im1 = limit(i-1, N)
jp1 = limit(j+1, N); jm1 = limit(j-1, N)
au = dfdα[II[i,j,1],II[i,j,1]] = (u[II[im1,j,1]] + u[II[ip1,j,1]] + u[II[i,jp1,1]] + u[II[i,jm1,1]] - 4u[II[i,j,1]])/dx^2
du[II[i,j,1]] = α[II[i,j,1]]*(au) + B[II[i,j,1]] + u[II[i,j,1]]^2*u[II[i,j,2]] - (A[II[i,j,1]] + 1)*u[II[i,j,1]] + brusselator_f(x, y, t)
end
for I in CartesianIndices((N, N))
i = I[1]
j = I[2]
ip1 = limit(i+1, N)
im1 = limit(i-1, N)
jp1 = limit(j+1, N)
jm1 = limit(j-1, N)
av = dfdα[II[i,j,2],II[i,j,2]] = (u[II[im1,j,2]] + u[II[ip1,j,2]] + u[II[i,jp1,2]] + u[II[i,jm1,2]] - 4u[II[i,j,2]])/dx^2
du[II[i,j,2]] = α[II[i,j,2]]*(av) + A[II[i,j,1]]*u[II[i,j,1]] - u[II[i,j,1]]^2*u[II[i,j,2]]
end
brusselator_jac(Jmat,u,p,t)
BLAS.gemm!('N', 'N', 1., Jmat, reshape(s, 2N*N, 4N*N), 1., dp)
copyto!(ds, vec(dp))
return nothing
end
end
end
u0 = init_brusselator_2d(xyd_brusselator)
p = [fill(3.4,N^2); fill(1.,N^2); fill(10.,2*N^2)]
brusselator_2d_loop, u0, p, brusselator_jac, ODEProblem{true, SciMLBase.FullSpecialize}(brusselator_comp, copy([u0;zeros((N^2*2)*(N^2*4))]), (0.,10.), p)
end
makebrusselator (generic function with 2 methods)
Differentiation Setups
function diffeq_sen(prob::DiffEqBase.DEProblem, args...; kwargs...)
diffeq_sen(prob.f, prob.u0, prob.tspan, prob.p, args...; kwargs...)
end
function auto_sen(prob::DiffEqBase.DEProblem, args...; kwargs...)
auto_sen(prob.f, prob.u0, prob.tspan, prob.p, args...; kwargs...)
end
function diffeq_sen(f, u0, tspan, p, alg=Tsit5(); sensalg=ForwardSensitivity(), kwargs...)
prob = ODEForwardSensitivityProblem(f,u0,tspan,p,sensalg)
sol = solve(prob,alg; save_everystep=false, kwargs...)
extract_local_sensitivities(sol, length(sol))[2]
end
function auto_sen(f, u0, tspan, p, alg=Tsit5(); kwargs...)
test_f(p) = begin
prob = ODEProblem{true, SciMLBase.FullSpecialize}(f,eltype(p).(u0),tspan,p)
solve(prob,alg; save_everystep=false, kwargs...)[end]
end
ForwardDiff.jacobian(test_f, p)
end
function numerical_sen(f,u0, tspan, p, alg=Tsit5(); kwargs...)
test_f(out,p) = begin
prob = ODEProblem{true, SciMLBase.FullSpecialize}(f,eltype(p).(u0),tspan,p)
copyto!(out, solve(prob,alg; kwargs...)[end])
end
J = Matrix{Float64}(undef,length(u0),length(p))
FiniteDiff.finite_difference_jacobian!(J, test_f, p, FiniteDiff.JacobianCache(p,Array{Float64}(undef,length(u0))))
return J
end
function diffeq_sen_l2(df, u0, tspan, p, t, alg=Tsit5();
abstol=1e-5, reltol=1e-7,
sensalg=InterpolatingAdjoint(), kwargs...)
prob = ODEProblem(df,u0,tspan,p)
sol = solve(prob, alg, sensealg=DiffEqBase.SensitivityADPassThrough(), abstol=abstol, reltol=reltol; kwargs...)
dg(out,u,p,t,i) = (out.=u.-1.0)
adjoint_sensitivities(sol,alg;t,abstol=abstol,dgdu_discrete = dg,
reltol=reltol,sensealg=sensalg)[2]
end
function auto_sen_l2(f, u0, tspan, p, t, alg=Tsit5(); diffalg=ReverseDiff.gradient, kwargs...)
test_f(p) = begin
prob = ODEProblem{true, SciMLBase.FullSpecialize}(f,eltype(p).(u0),tspan,p)
sol = solve(prob,alg; sensealg=DiffEqBase.SensitivityADPassThrough(), kwargs...)(t)
sum(sol.u) do x
sum(z->(1-z)^2/2, x)
end
end
diffalg(test_f, p)
end
function numerical_sen_l2(f, u0, tspan, p, t, alg=Tsit5(); kwargs...)
test_f(p) = begin
prob = ODEProblem(f,eltype(p).(u0),tspan,p)
sol = solve(prob,alg; kwargs...)(t)
sum(sol.u) do x
sum(z->(1-z)^2/2, x)
end
end
FiniteDiff.finite_difference_gradient(test_f, p, Val{:central})
end
numerical_sen_l2 (generic function with 2 methods)
_adjoint_methods = ntuple(3) do ii
Alg = (InterpolatingAdjoint, QuadratureAdjoint, BacksolveAdjoint)[ii]
(
user = Alg(autodiff=false,autojacvec=false), # user Jacobian
adjc = Alg(autodiff=true,autojacvec=false), # AD Jacobian
advj = Alg(autodiff=true,autojacvec=EnzymeVJP()), # AD vJ
)
end |> NamedTuple{(:interp, :quad, :backsol)}
@isdefined(ADJOINT_METHODS) || (const ADJOINT_METHODS = mapreduce(collect, vcat, _adjoint_methods))
9-element Vector{SciMLBase.AbstractAdjointSensitivityAlgorithm{0, AD, Val{:
central}} where AD}:
SciMLSensitivity.InterpolatingAdjoint{0, false, Val{:central}, Bool}(false
, false, false)
SciMLSensitivity.InterpolatingAdjoint{0, true, Val{:central}, Bool}(false,
false, false)
SciMLSensitivity.InterpolatingAdjoint{0, true, Val{:central}, SciMLSensiti
vity.EnzymeVJP}(SciMLSensitivity.EnzymeVJP(0), false, false)
SciMLSensitivity.QuadratureAdjoint{0, false, Val{:central}, Bool}(false, 1
.0e-6, 0.001)
SciMLSensitivity.QuadratureAdjoint{0, true, Val{:central}, Bool}(false, 1.
0e-6, 0.001)
SciMLSensitivity.QuadratureAdjoint{0, true, Val{:central}, SciMLSensitivit
y.EnzymeVJP}(SciMLSensitivity.EnzymeVJP(0), 1.0e-6, 0.001)
SciMLSensitivity.BacksolveAdjoint{0, false, Val{:central}, Bool}(false, tr
ue, false)
SciMLSensitivity.BacksolveAdjoint{0, true, Val{:central}, Bool}(false, tru
e, false)
SciMLSensitivity.BacksolveAdjoint{0, true, Val{:central}, SciMLSensitivity
.EnzymeVJP}(SciMLSensitivity.EnzymeVJP(0), true, false)
Run Forward Mode Benchmarks
These are testing for the construction of the full Jacobian.
forward_lv = let
@info "Running the Lotka-Volterra model:"
@info " Running compile-time CSA"
t1 = @belapsed solve($lvprob, $(Tsit5()); $tols...)
@info " Running DSA"
t2 = @belapsed auto_sen($lvdf, $u0, $tspan, $p, $(Tsit5()); $tols...)
@info " Running CSA user-Jacobian"
t3 = @belapsed diffeq_sen($lvdf_with_jacobian, $u0, $tspan, $p, $(Tsit5()); sensalg=ForwardSensitivity(autodiff=false, autojacvec=false), $tols...)
@info " Running AD-Jacobian"
t4 = @belapsed diffeq_sen($lvdf, $u0, $tspan, $p, $(Tsit5()); sensalg=ForwardSensitivity(autojacvec=false), $tols...)
@info " Running AD-Jv seeding"
t5 = @belapsed diffeq_sen($lvdf, $u0, $tspan, $p, $(Tsit5()); sensalg=ForwardSensitivity(autojacvec=true), $tols...)
@info " Running numerical differentiation"
t6 = @belapsed numerical_sen($lvdf, $u0, $tspan, $p, $(Tsit5()); $tols...)
print('\n')
[t1, t2, t3, t4, t5, t6]
end
6-element Vector{Float64}:
7.585e-5
3.7719e-5
0.000311147
0.000423807
0.000287728
0.000184449
forward_bruss = let
@info "Running the Brusselator model:"
n = 5
# Run low tolerance to test correctness
bfun, b_u0, b_p, brusselator_jac, brusselator_comp = makebrusselator(n)
sol1 = @time numerical_sen(bfun, b_u0, (0.,10.), b_p, Rodas5(), abstol=1e-5,reltol=1e-7);
sol2 = @time auto_sen(bfun, b_u0, (0.,10.), b_p, Rodas5(), abstol=1e-5,reltol=1e-7);
@test sol1 ≈ sol2 atol=1e-2
sol3 = @time diffeq_sen(bfun, b_u0, (0.,10.), b_p, Rodas5(autodiff=false), abstol=1e-5,reltol=1e-7);
@test sol1 ≈ hcat(sol3...) atol=1e-3
sol4 = @time diffeq_sen(ODEFunction{true, SciMLBase.FullSpecialize}(bfun, jac=brusselator_jac), b_u0, (0.,10.), b_p, Rodas5(autodiff=false), abstol=1e-5,reltol=1e-7, sensalg=ForwardSensitivity(autodiff=false, autojacvec=false));
@test sol1 ≈ hcat(sol4...) atol=1e-2
sol5 = @time solve(brusselator_comp, Rodas5(autodiff=false), abstol=1e-5,reltol=1e-7,);
@test sol1 ≈ reshape(sol5[end][2n*n+1:end], 2n*n, 4n*n) atol=1e-3
# High tolerance to benchmark
@info " Running compile-time CSA"
t1 = @belapsed solve($brusselator_comp, $(Rodas5(autodiff=false)); $tols...);
@info " Running DSA"
t2 = @belapsed auto_sen($bfun, $b_u0, $((0.,10.)), $b_p, $(Rodas5()); $tols...);
@info " Running CSA user-Jacobian"
t3 = @belapsed diffeq_sen($(ODEFunction{true, SciMLBase.FullSpecialize}(bfun, jac=brusselator_jac)), $b_u0, $((0.,10.)), $b_p, $(Rodas5(autodiff=false)); sensalg=ForwardSensitivity(autodiff=false, autojacvec=false), $tols...);
@info " Running AD-Jacobian"
t4 = @belapsed diffeq_sen($bfun, $b_u0, $((0.,10.)), $b_p, $(Rodas5(autodiff=false)); sensalg=ForwardSensitivity(autojacvec=false), $tols...);
@info " Running AD-Jv seeding"
t5 = @belapsed diffeq_sen($bfun, $b_u0, $((0.,10.)), $b_p, $(Rodas5(autodiff=false)); sensalg=ForwardSensitivity(autojacvec=true), $tols...);
@info " Running numerical differentiation"
t6 = @belapsed numerical_sen($bfun, $b_u0, $((0.,10.)), $b_p, $(Rodas5()); $tols...);
print('\n')
[t1, t2, t3, t4, t5, t6]
end
5.260279 seconds (10.06 M allocations: 691.162 MiB, 5.42% gc time, 92.74%
compilation time)
16.438322 seconds (19.12 M allocations: 1.436 GiB, 6.40% gc time, 86.01% c
ompilation time)
76.284152 seconds (5.84 M allocations: 765.213 MiB, 0.24% gc time, 6.69% c
ompilation time)
103.874826 seconds (7.30 M allocations: 1.015 GiB, 2.62% gc time, 4.96% com
pilation time)
62.315488 seconds (4.23 M allocations: 816.423 MiB, 1.77% gc time, 3.82% c
ompilation time)
6-element Vector{Float64}:
52.110420308
1.753122303
87.810682803
64.591683244
72.974257057
0.341829561
forward_pollution = let
@info "Running the pollution model:"
pcomp, pu0, pp, pcompu0 = make_pollution()
ptspan = (0.,60.)
@info " Running compile-time CSA"
t1 = 0#@belapsed solve($(ODEProblem(pcomp, pcompu0, ptspan, pp)), $(Rodas5(autodiff=false)),);
@info " Running DSA"
t2 = @belapsed auto_sen($(ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), $pu0, $ptspan, $pp, $(Rodas5()); $tols...);
@info " Running CSA user-Jacobian"
t3 = @belapsed diffeq_sen($(ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f, jac=pollution.jac)), $pu0, $ptspan, $pp, $(Rodas5(autodiff=false)); sensalg=ForwardSensitivity(autodiff=false, autojacvec=false), $tols...);
@info " Running AD-Jacobian"
t4 = @belapsed diffeq_sen($(ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), $pu0, $ptspan, $pp, $(Rodas5(autodiff=false)); sensalg=ForwardSensitivity(autojacvec=false), $tols...);
@info " Running AD-Jv seeding"
t5 = @belapsed diffeq_sen($(ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), $pu0, $ptspan, $pp, $(Rodas5(autodiff=false)); sensalg=ForwardSensitivity(autojacvec=true), $tols...);
@info " Running numerical differentiation"
t6 = @belapsed numerical_sen($(ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), $pu0, $ptspan, $pp, $(Rodas5()); $tols...);
print('\n')
[t1, t2, t3, t4, t5, t6]
end
6-element Vector{Float64}:
0.0
0.007768163
0.435246767
0.432378448
0.421477738
0.006342933
forward_pkpd = let
@info "Running the PKPD model:"
#sol1 = solve(pkpdcompprob, Tsit5(),abstol=1e-5,reltol=1e-7,callback=pkpdcb,tstops=0:24:240,)[end][6:end]
sol2 = vec(auto_sen(pkpdprob, Tsit5(),abstol=1e-5,reltol=1e-7,callback=pkpdcb,tstops=0:24:240))
sol3 = vec(hcat(diffeq_sen(pkpdprob, Tsit5(),abstol=1e-5,reltol=1e-7,callback=pkpdcb,tstops=0:24:240)...))
#@test sol1 ≈ sol2 atol=1e-3
@test sol2 ≈ sol3 atol=1e-3
@info " Running compile-time CSA"
#t1 = @belapsed solve($pkpdcompprob, $(Tsit5()),callback=$pkpdcb,tstops=0:24:240,);
@info " Running DSA"
t2 = @belapsed auto_sen($(pkpdf.f), $pkpdu0, $pkpdtspan, $pkpdp, $(Tsit5()); callback=$pkpdcb,tstops=0:24:240, $tols...);
@info " Running CSA user-Jacobian"
t3 = @belapsed diffeq_sen($(ODEFunction{true, SciMLBase.FullSpecialize}(pkpdf.f, jac=pkpdf.jac)), $pkpdu0, $pkpdtspan, $pkpdp, $(Tsit5());callback=$pkpdcb,tstops=0:24:240, sensalg=ForwardSensitivity(autodiff=false, autojacvec=false), $tols...);
@info " Running AD-Jacobian"
t4 = @belapsed diffeq_sen($(pkpdf.f), $pkpdu0, $pkpdtspan, $pkpdp, $(Tsit5()); callback=$pkpdcb,tstops=0:24:240,
sensalg=ForwardSensitivity(autojacvec=false), $tols...);
@info " Running AD-Jv seeding"
t5 = @belapsed diffeq_sen($(pkpdf.f), $pkpdu0, $pkpdtspan, $pkpdp, $(Tsit5());callback=$pkpdcb,tstops=0:24:240,
sensalg=ForwardSensitivity(autojacvec=true), $tols...);
@info " Running numerical differentiation"
t6 = @belapsed numerical_sen($(pkpdf.f), $pkpdu0, $pkpdtspan, $pkpdp, $(Tsit5()); callback=$pkpdcb,tstops=0:24:240, $tols...);
print('\n')
[0, t2, t3, t4, t5, t6]
end
6-element Vector{Float64}:
0.0
0.001661177
0.00669911
0.005292711
0.007619544
0.006508272
forward_methods = ["Compile-time CSA", "DSA", "CSA user-Jacobian", "AD-Jacobian", "AD-Jv seeding", "Numerical Differentiation"]
forward_timings = DataFrame(methods=forward_methods, LV=forward_lv, Bruss=forward_bruss, Pollution=forward_pollution, PKPD=forward_pkpd)
display(forward_timings)
6×5 DataFrame
Row │ methods LV Bruss Pollution PKPD
⋯
│ String Float64 Float64 Float64 Float6
4 ⋯
─────┼─────────────────────────────────────────────────────────────────────
─────
1 │ Compile-time CSA 7.585e-5 52.1104 0.0 0.0
⋯
2 │ DSA 3.7719e-5 1.75312 0.00776816 0.0016
611
3 │ CSA user-Jacobian 0.000311147 87.8107 0.435247 0.0066
991
4 │ AD-Jacobian 0.000423807 64.5917 0.432378 0.0052
927
5 │ AD-Jv seeding 0.000287728 72.9743 0.421478 0.0076
195 ⋯
6 │ Numerical Differentiation 0.000184449 0.34183 0.00634293 0.0065
082
1 column om
itted
Run Adjoint Benchmarks
Adjoint requires a slightly different setup even with forward mode ADs since it requires a loss function choice. For that we simply take the L2 norm of the solution.
adjoint_lv = let
@info "Running the Lotka-Volerra model:"
lvu0 = [1.,1.]; lvtspan = (0.0, 10.0); lvp = [1.5,1.0,3.0];
lvt = 0:0.5:10
@time lsol1 = auto_sen_l2(lvdf, lvu0, lvtspan, lvp, lvt, (Tsit5()); diffalg=(ForwardDiff.gradient), tols...);
@time lsol2 = auto_sen_l2(lvdf, lvu0, lvtspan, lvp, lvt, (Tsit5()); diffalg=(ReverseDiff.gradient), tols...);
@time lsol3 = map(ADJOINT_METHODS) do alg
f = SciMLSensitivity.alg_autodiff(alg) ? lvdf : lvdf_with_jacobian
diffeq_sen_l2(f, lvu0, lvtspan, lvp, lvt, (Tsit5()); sensalg=alg, tols...)
end
@time lsol4 = numerical_sen_l2(lvdf, lvu0, lvtspan, lvp, lvt, Tsit5(); tols...);
@test maximum(abs, lsol1 .- lsol2)/maximum(abs, lsol1) < 0.2
@test all(i -> maximum(abs, lsol1 .- lsol3[i]')/maximum(abs, lsol1) < 0.2, eachindex(ADJOINT_METHODS))
@test maximum(abs, lsol1 .- lsol4)/maximum(abs, lsol1) < 0.2
t1 = @belapsed auto_sen_l2($lvdf, $lvu0, $lvtspan, $lvp, $lvt, $(Tsit5()); diffalg=$(ForwardDiff.gradient), $tols...);
t2 = @belapsed auto_sen_l2($lvdf, $lvu0, $lvtspan, $lvp, $lvt, $(Tsit5()); diffalg=$(ReverseDiff.gradient), $tols...);
t3 = map(ADJOINT_METHODS) do alg
f = SciMLSensitivity.alg_autodiff(alg) ? lvdf : lvdf_with_jacobian
@belapsed diffeq_sen_l2($f, $lvu0, $lvtspan, $lvp, $lvt, $(Tsit5()); sensalg=$alg, $tols...);
end
t4 = @belapsed numerical_sen_l2($lvdf, $lvu0, $lvtspan, $lvp, $lvt, $(Tsit5()); $tols...);
[t1; t2; t3; t4]
end
2.665906 seconds (5.96 M allocations: 401.429 MiB, 4.09% gc time, 99.94%
compilation time)
4.101357 seconds (6.85 M allocations: 467.553 MiB, 3.77% gc time, 99.82%
compilation time)
53.767595 seconds (82.43 M allocations: 5.496 GiB, 3.39% gc time, 99.79% c
ompilation time: <1% of which was recompilation)
0.507162 seconds (615.78 k allocations: 41.886 MiB, 99.56% compilation ti
me)
12-element Vector{Float64}:
9.8019e-5
0.004172429
0.000444077
0.000545935
0.000330187
0.001077622
0.001199211
0.001709738
0.000756614
0.000920413
0.000493336
0.000604696
adjoint_bruss = let
@info "Running the Brusselator model:"
bt = 0:0.1:10
tspan = (0.0, 10.0)
n = 5
bfun, b_u0, b_p, brusselator_jac, brusselator_comp = makebrusselator(n)
@time bsol1 = auto_sen_l2(bfun, b_u0, tspan, b_p, bt, (Rodas5()); diffalg=(ForwardDiff.gradient), tols...);
#@time bsol2 = auto_sen_l2(bfun, b_u0, tspan, b_p, bt, (Rodas5(autodiff=false)); diffalg=(ReverseDiff.gradient), tols...);
#@test maximum(abs, bsol1 .- bsol2)/maximum(abs, bsol1) < 1e-2
@time bsol3 = map(ADJOINT_METHODS) do alg
@info "Running $alg"
f = SciMLSensitivity.alg_autodiff(alg) ? bfun : ODEFunction{true, SciMLBase.FullSpecialize}(bfun, jac=brusselator_jac)
solver = Rodas5(autodiff=false)
diffeq_sen_l2(f, b_u0, tspan, b_p, bt, solver, reltol=1e-7; sensalg=alg, tols...)
end
@time bsol4 = numerical_sen_l2(bfun, b_u0, tspan, b_p, bt, (Rodas5()); tols...);
# NOTE: backsolve gives unstable results!!!
@test all(i->maximum(abs, bsol1 .- bsol3[i]')/maximum(abs, bsol1) < 4e-2, eachindex(ADJOINT_METHODS)[1:2end÷3])
@test all(i->maximum(abs, bsol1 .- bsol3[i]')/maximum(abs, bsol1) >= 4e-2, eachindex(ADJOINT_METHODS)[2end÷3+1:end])
@test maximum(abs, bsol1 .- bsol4)/maximum(abs, bsol1) < 2e-2
t1 = @belapsed auto_sen_l2($bfun, $b_u0, $tspan, $b_p, $bt, $(Rodas5()); diffalg=$(ForwardDiff.gradient), $tols...);
#t2 = @belapsed auto_sen_l2($bfun, $b_u0, $tspan, $b_p, $bt, $(Rodas5(autodiff=false)); diffalg=$(ReverseDiff.gradient), $tols...);
t2 = NaN
t3 = map(ADJOINT_METHODS[1:2end÷3]) do alg
@info "Running $alg"
f = SciMLSensitivity.alg_autodiff(alg) ? bfun : ODEFunction{true, SciMLBase.FullSpecialize}(bfun, jac=brusselator_jac)
solver = Rodas5(autodiff=false)
@elapsed diffeq_sen_l2(f, b_u0, tspan, b_p, bt, solver; sensalg=alg, tols...);
end
t3 = [t3; fill(NaN, length(ADJOINT_METHODS)÷3)]
t4 = @belapsed numerical_sen_l2($bfun, $b_u0, $tspan, $b_p, $bt, $(Rodas5()); $tols...);
[t1; t2; t3; t4]
end
15.664942 seconds (17.92 M allocations: 1.374 GiB, 2.94% gc time, 87.87% c
ompilation time)
75.557440 seconds (69.08 M allocations: 4.416 GiB, 1.60% gc time, 88.18% c
ompilation time)
3.306743 seconds (2.75 M allocations: 226.488 MiB, 1.68% gc time, 70.18%
compilation time)
12-element Vector{Float64}:
1.762203879
NaN
2.865565296
1.455263483
0.099712688
0.170598081
0.261438242
0.056026099
NaN
NaN
NaN
0.937964967
adjoint_pollution = let
@info "Running the Pollution model:"
pcomp, pu0, pp, pcompu0 = make_pollution();
ptspan = (0.0, 60.0)
pts = 0:0.5:60
@time psol1 = auto_sen_l2((ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), pu0, ptspan, pp, pts, (Rodas5(autodiff=false)); diffalg=(ForwardDiff.gradient), tols...);
#@time psol2 = auto_sen_l2((ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), pu0, ptspan, pp, pts, (Rodas5(autodiff=false)); diffalg=(ReverseDiff.gradient), tols...);
#@test maximum(abs, psol1 .- psol2)/maximum(abs, psol1) < 1e-2
@time psol3 = map(ADJOINT_METHODS) do alg
@info "Running $alg"
f = SciMLSensitivity.alg_autodiff(alg) ? pollution.f : ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f, jac=pollution.jac)
solver = Rodas5(autodiff=false)
diffeq_sen_l2(f, pu0, ptspan, pp, pts, solver; sensalg=alg, tols...);
end
@time psol4 = numerical_sen_l2((ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), pu0, ptspan, pp, pts, (Rodas5(autodiff=false)); tols...);
# NOTE: backsolve gives unstable results!!!
@test all(i->maximum(abs, psol1 .- psol3[i]')/maximum(abs, psol1) < 1e-2, eachindex(ADJOINT_METHODS)[1:2end÷3])
@test all(i->maximum(abs, psol1 .- psol3[i]')/maximum(abs, psol1) >= 1e-2, eachindex(ADJOINT_METHODS)[2end÷3+1:end])
@test maximum(abs, psol1 .- psol4)/maximum(abs, psol1) < 1e-2
t1 = @belapsed auto_sen_l2($(ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), $pu0, $ptspan, $pp, $pts, $(Rodas5(autodiff=false)); diffalg=$(ForwardDiff.gradient), $tols...);
#t2 = @belapsed auto_sen_l2($(ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), $pu0, $ptspan, $pp, $pts, $(Rodas5(autodiff=false)); diffalg=$(ReverseDiff.gradient), $tols...);
t2 = NaN
t3 = map(ADJOINT_METHODS[1:2end÷3]) do alg
@info "Running $alg"
f = SciMLSensitivity.alg_autodiff(alg) ? pollution.f : ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f, jac=pollution.jac)
solver = Rodas5(autodiff=false)
@elapsed diffeq_sen_l2(f, pu0, ptspan, pp, pts, solver; sensalg=alg, tols...);
end
t3 = [t3; fill(NaN, length(ADJOINT_METHODS)÷3)]
t4 = @belapsed numerical_sen_l2($(ODEFunction{true, SciMLBase.FullSpecialize}(pollution.f)), $pu0, $ptspan, $pp, $pts, $(Rodas5(autodiff=false)); $tols...);
[t1; t2; t3; t4]
end
15.105752 seconds (17.83 M allocations: 1.110 GiB, 2.73% gc time, 99.95% c
ompilation time)
62.069172 seconds (68.81 M allocations: 3.899 GiB, 1.96% gc time, 91.43% c
ompilation time)
2.324609 seconds (2.55 M allocations: 156.528 MiB, 2.81% gc time, 99.31%
compilation time)
12-element Vector{Float64}:
0.005063724
NaN
0.823905908
1.288597705
0.194944666
0.19574358
0.573064572
0.3524813
NaN
NaN
NaN
0.014632476
adjoint_pkpd = let
@info "Running the PKPD model:"
pts = 0:0.5:50
# need to use lower tolerances to avoid running into the complex domain because of exponentiation
pkpdsol1 = @time auto_sen_l2((pkpdf.f), pkpdu0, pkpdtspan, pkpdp, pts, (Tsit5()); callback=pkpdcb, tstops=0:24:240,
diffalg=(ForwardDiff.gradient), tols...);
pkpdsol2 = @time auto_sen_l2((pkpdf.f), pkpdu0, pkpdtspan, pkpdp, pts, (Tsit5()); callback=pkpdcb, tstops=0:24:240,
diffalg=(ReverseDiff.gradient), tols...);
pkpdsol3 = @time map(ADJOINT_METHODS[1:2end÷3]) do alg
f = SciMLSensitivity.alg_autodiff(alg) ? pkpdf.f : ODEFunction{true, SciMLBase.FullSpecialize}(pkpdf.f, jac=pkpdf.jac)
diffeq_sen_l2(f, pkpdu0, pkpdtspan, pkpdp, pts, (Tsit5()); sensalg=alg,
callback=pkpdcb, tstops=0:24:240, tols...);
end
pkpdsol4 = @time numerical_sen_l2((ODEFunction{true, SciMLBase.FullSpecialize}(pkpdf.f)), pkpdu0, pkpdtspan, pkpdp, pts, (Tsit5());
callback=pkpdcb, tstops=0:24:240, tols...);
@test maximum(abs, pkpdsol1 .- pkpdsol2)/maximum(abs, pkpdsol1) < 0.2
@test all(i->maximum(abs, pkpdsol1 .- pkpdsol3[i]')/maximum(abs, pkpdsol1) < 0.2, eachindex(ADJOINT_METHODS)[1:2end÷3])
@test maximum(abs, pkpdsol1 .- pkpdsol4)/maximum(abs, pkpdsol1) < 0.2
t1 = @belapsed auto_sen_l2($(pkpdf.f), $pkpdu0, $pkpdtspan, $pkpdp, $pts, $(Tsit5()); callback=pkpdcb, tstops=0:24:240,
diffalg=$(ForwardDiff.gradient), $tols...);
t2 = @belapsed auto_sen_l2($(pkpdf.f), $pkpdu0, $pkpdtspan, $pkpdp, $pts, $(Tsit5()); callback=pkpdcb, tstops=0:24:240,
diffalg=$(ReverseDiff.gradient), $tols...);
t3 = map(ADJOINT_METHODS[1:2end÷3]) do alg
f = SciMLSensitivity.alg_autodiff(alg) ? pkpdf.f : ODEFunction{true, SciMLBase.FullSpecialize}(pkpdf.f, jac=pkpdf.jac)
@belapsed diffeq_sen_l2($f, $pkpdu0, $pkpdtspan, $pkpdp, $pts, $(Tsit5()); tstops=0:24:240,
callback=pkpdcb, sensalg=$alg, tols...);
end
t3 = [t3; fill(NaN, length(ADJOINT_METHODS)÷3)]
t4 = @belapsed numerical_sen_l2($(ODEFunction{true, SciMLBase.FullSpecialize}(pkpdf.f)), $pkpdu0, $pkpdtspan, $pkpdp, $pts, $(Tsit5()); tstops=0:24:240,
callback=$pkpdcb, $tols...);
[t1; t2; t3; t4]
end
2.981166 seconds (5.70 M allocations: 385.289 MiB, 2.68% gc time, 99.86%
compilation time)
3.108587 seconds (5.63 M allocations: 315.093 MiB, 4.63% gc time, 93.17%
compilation time)
26.707226 seconds (41.72 M allocations: 2.737 GiB, 2.91% gc time, 99.76% c
ompilation time)
0.283416 seconds (441.84 k allocations: 33.775 MiB, 15.45% gc time, 77.89
% compilation time)
12-element Vector{Float64}:
0.002349303
0.119229007
0.006465284
0.005010835
0.001725568
0.006247995
0.005750698
0.006056756
NaN
NaN
NaN
0.015676778
adjoint_methods = ["ForwardDiff", "ReverseDiff",
"InterpolatingAdjoint User Jac", "InterpolatingAdjoint AD Jac", "InterpolatingAdjoint v'J",
"QuadratureAdjoint User Jac", "QuadratureAdjoint AD Jac", "QuadratureAdjoint v'J",
"BacksolveAdjoint User Jac", "BacksolveAdjoint AD Jac", "BacksolveAdjoint v'J",
"Numerical Differentiation"]
adjoint_timings = DataFrame(methods=adjoint_methods, LV=adjoint_lv, Bruss=adjoint_bruss, Pollution=adjoint_pollution, PKPD=adjoint_pkpd)
Markdown.parse(PrettyTables.pretty_table(String, adjoint_timings; backend=Val(:markdown), header=names(adjoint_timings)))
methods | LV | Bruss | Pollution | PKPD |
---|---|---|---|---|
ForwardDiff | 9.8019e-5 | 1.7622 | 0.00506372 | 0.0023493 |
ReverseDiff | 0.00417243 | NaN | NaN | 0.119229 |
InterpolatingAdjoint User Jac | 0.000444077 | 2.86557 | 0.823906 | 0.00646528 |
InterpolatingAdjoint AD Jac | 0.000545935 | 1.45526 | 1.2886 | 0.00501083 |
InterpolatingAdjoint v'J | 0.000330187 | 0.0997127 | 0.194945 | 0.00172557 |
QuadratureAdjoint User Jac | 0.00107762 | 0.170598 | 0.195744 | 0.006248 |
QuadratureAdjoint AD Jac | 0.00119921 | 0.261438 | 0.573065 | 0.0057507 |
QuadratureAdjoint v'J | 0.00170974 | 0.0560261 | 0.352481 | 0.00605676 |
BacksolveAdjoint User Jac | 0.000756614 | NaN | NaN | NaN |
BacksolveAdjoint AD Jac | 0.000920413 | NaN | NaN | NaN |
BacksolveAdjoint v'J | 0.000493336 | NaN | NaN | NaN |
Numerical Differentiation | 0.000604696 | 0.937965 | 0.0146325 | 0.0156768 |
Appendix
Appendix
These benchmarks are a part of the SciMLBenchmarks.jl repository, found at: https://github.com/SciML/SciMLBenchmarks.jl. For more information on high-performance scientific machine learning, check out the SciML Open Source Software Organization https://sciml.ai.
To locally run this benchmark, do the following commands:
using SciMLBenchmarks
SciMLBenchmarks.weave_file("benchmarks/AutomaticDifferentiation","SimpleODEAD.jmd")
Computer Information:
Julia Version 1.10.9
Commit 5595d20a287 (2025-03-10 12:51 UTC)
Build Info:
Official https://julialang.org/ release
Platform Info:
OS: Linux (x86_64-linux-gnu)
CPU: 128 × AMD EPYC 7502 32-Core Processor
WORD_SIZE: 64
LIBM: libopenlibm
LLVM: libLLVM-15.0.7 (ORCJIT, znver2)
Threads: 1 default, 0 interactive, 1 GC (on 128 virtual cores)
Environment:
JULIA_CPU_THREADS = 128
JULIA_DEPOT_PATH = /cache/julia-buildkite-plugin/depots/5b300254-1738-4989-ae0a-f4d2d937f953
Package Information:
Status `/cache/build/exclusive-amdci1-0/julialang/scimlbenchmarks-dot-jl/benchmarks/AutomaticDifferentiation/Project.toml`
⌃ [6e4b80f9] BenchmarkTools v1.5.0
⌃ [a93c6f00] DataFrames v1.6.1
⌃ [1313f7d8] DataFramesMeta v0.15.3
⌅ [a0c0ee7d] DifferentiationInterface v0.5.9
⌅ [a82114a7] DifferentiationInterfaceTest v0.5.0
⌅ [7da242da] Enzyme v0.12.25
⌃ [6a86dc24] FiniteDiff v2.23.1
⌅ [f6369f11] ForwardDiff v0.10.36
⌃ [1dea7af3] OrdinaryDiffEq v6.86.0
⌃ [65888b18] ParameterizedFunctions v5.17.0
⌃ [91a5bcdd] Plots v1.40.5
⌃ [08abe8d2] PrettyTables v2.3.2
[37e2e3b7] ReverseDiff v1.15.3
[31c91b34] SciMLBenchmarks v0.1.3
⌃ [1ed8b502] SciMLSensitivity v7.64.0
⌃ [90137ffa] StaticArrays v1.9.7
⌃ [07d77754] Tapir v0.2.26
⌃ [9f7883ad] Tracker v0.2.34
⌅ [e88e6eb3] Zygote v0.6.70
[37e2e46d] LinearAlgebra
[d6f4376e] Markdown
[de0858da] Printf
[8dfed614] Test
Info Packages marked with ⌃ and ⌅ have new versions available. Those with ⌃ may be upgradable, but those with ⌅ are restricted by compatibility constraints from upgrading. To see why use `status --outdated`
Warning The project dependencies or compat requirements have changed since the manifest was last resolved. It is recommended to `Pkg.resolve()` or consider `Pkg.update()` if necessary.
And the full manifest:
Status `/cache/build/exclusive-amdci1-0/julialang/scimlbenchmarks-dot-jl/benchmarks/AutomaticDifferentiation/Manifest.toml`
⌃ [47edcb42] ADTypes v1.6.1
[621f4979] AbstractFFTs v1.5.0
[1520ce14] AbstractTrees v0.4.5
⌃ [7d9f7c33] Accessors v0.1.37
⌃ [79e6a3ab] Adapt v4.0.4
[66dad0bd] AliasTables v1.1.3
[ec485272] ArnoldiMethod v0.4.0
⌃ [4fba245c] ArrayInterface v7.12.0
⌃ [4c555306] ArrayLayouts v1.10.2
⌅ [a9b6321e] Atomix v0.1.0
⌃ [6e4b80f9] BenchmarkTools v1.5.0
⌃ [e2ed5e7c] Bijections v0.1.7
[d1d4a3ce] BitFlags v0.1.9
[62783981] BitTwiddlingConvenienceFunctions v0.1.6
[fa961155] CEnum v0.5.0
[2a0fbf3d] CPUSummary v0.2.6
[00ebfdb7] CSTParser v3.4.3
⌃ [49dc2e85] Calculus v0.5.1
⌃ [7057c7e9] Cassette v0.3.13
[8be319e6] Chain v0.6.0
⌃ [082447d4] ChainRules v1.69.0
⌃ [d360d2e6] ChainRulesCore v1.24.0
⌃ [0ca39b1e] Chairmarks v1.2.1
[fb6a15b2] CloseOpenIntervals v0.1.13
⌃ [da1fd8a2] CodeTracking v1.3.5
⌃ [944b1d66] CodecZlib v0.7.5
⌃ [35d6a980] ColorSchemes v3.26.0
⌅ [3da002f7] ColorTypes v0.11.5
⌅ [c3611d14] ColorVectorSpace v0.10.0
⌅ [5ae59095] Colors v0.12.11
[861a8166] Combinatorics v1.0.2
⌅ [a80b9123] CommonMark v0.8.12
[38540f10] CommonSolve v0.2.4
⌃ [bbf7d656] CommonSubexpressions v0.3.0
[f70d9fcc] CommonWorldInvalidations v1.0.0
⌃ [34da2185] Compat v4.15.0
⌃ [b0b7db55] ComponentArrays v0.15.14
[b152e2b5] CompositeTypes v0.1.4
[a33af91c] CompositionsBase v0.1.2
[2569d6c7] ConcreteStructs v0.2.3
⌃ [f0e56b4a] ConcurrentUtilities v2.4.2
[8f4d0f93] Conda v1.10.2
⌅ [187b0558] ConstructionBase v1.5.6
[d38c429a] Contour v0.6.3
[adafc99b] CpuId v0.3.1
[a8cc5b0e] Crayons v4.1.1
[9a962f9c] DataAPI v1.16.0
⌃ [a93c6f00] DataFrames v1.6.1
⌃ [1313f7d8] DataFramesMeta v0.15.3
⌃ [864edb3b] DataStructures v0.18.20
[e2d170a0] DataValueInterfaces v1.0.0
[8bb1440f] DelimitedFiles v1.9.1
⌃ [2b5f629d] DiffEqBase v6.151.5
⌅ [459566f4] DiffEqCallbacks v3.6.2
⌃ [77a26b50] DiffEqNoiseProcess v5.22.0
[163ba53b] DiffResults v1.1.0
[b552c78f] DiffRules v1.15.1
[de460e47] DiffTests v0.1.2
⌅ [a0c0ee7d] DifferentiationInterface v0.5.9
⌅ [a82114a7] DifferentiationInterfaceTest v0.5.0
⌃ [b4f34e82] Distances v0.10.11
⌃ [31c24e10] Distributions v0.25.109
[ffbed154] DocStringExtensions v0.9.3
⌃ [5b8099bc] DomainSets v0.7.14
⌃ [fa6b7ba4] DualNumbers v0.6.8
⌅ [7c1d4256] DynamicPolynomials v0.5.7
⌅ [06fc5a27] DynamicQuantities v0.13.2
[da5c29d0] EllipsisNotation v1.8.0
[4e289a0a] EnumX v1.0.4
⌅ [7da242da] Enzyme v0.12.25
⌅ [f151be2c] EnzymeCore v0.7.7
⌃ [460bff9d] ExceptionUnwrapping v0.1.10
⌃ [d4d017d3] ExponentialUtilities v1.26.1
[e2ba6199] ExprTools v0.1.10
⌃ [c87230d0] FFMPEG v0.4.1
⌃ [7034ab61] FastBroadcast v0.3.4
[9aa1b823] FastClosures v0.3.2
[29a986be] FastLapackInterface v2.0.4
⌃ [1a297f60] FillArrays v1.11.0
⌃ [64ca27bc] FindFirstFunctions v1.2.0
⌃ [6a86dc24] FiniteDiff v2.23.1
[53c48c17] FixedPointNumbers v0.8.5
[1fa38f19] Format v1.3.7
⌅ [f6369f11] ForwardDiff v0.10.36
[f62d2435] FunctionProperties v0.1.2
[069b7b12] FunctionWrappers v1.1.3
[77dc65aa] FunctionWrappersWrappers v0.1.3
⌅ [d9f16b24] Functors v0.4.11
⌅ [0c68f7d7] GPUArrays v10.3.0
⌅ [46192b85] GPUArraysCore v0.1.6
⌅ [61eb1bfa] GPUCompiler v0.26.7
⌃ [28b8d3ca] GR v0.73.7
[c145ed77] GenericSchur v0.5.4
[d7ba0133] Git v1.3.1
[c27321d9] Glob v1.3.1
⌃ [86223c79] Graphs v1.11.2
[42e2da0e] Grisu v1.0.2
⌃ [cd3eb016] HTTP v1.10.8
[eafb193a] Highlights v0.5.3
[3e5b6fbb] HostCPUFeatures v0.1.17
⌃ [34004b35] HypergeometricFunctions v0.3.23
⌃ [7073ff75] IJulia v1.25.0
[7869d1d1] IRTools v0.4.14
[615f187c] IfElse v0.1.1
[d25df0c9] Inflate v0.1.5
⌃ [842dd82b] InlineStrings v1.4.2
[8197267c] IntervalSets v0.7.10
⌃ [3587e190] InverseFunctions v0.1.15
⌃ [41ab1584] InvertedIndices v1.3.0
⌃ [92d709cd] IrrationalConstants v0.2.2
[82899510] IteratorInterfaceExtensions v1.0.0
⌃ [c3a54625] JET v0.9.6
⌅ [27aeb0d3] JLArrays v0.1.5
⌃ [1019f520] JLFzf v0.1.7
⌃ [692b3bcd] JLLWrappers v1.5.0
[682c06a0] JSON v0.21.4
⌃ [98e50ef6] JuliaFormatter v1.0.58
⌃ [aa1ae85d] JuliaInterpreter v0.9.32
⌃ [ccbc3e58] JumpProcesses v9.11.1
[ef3ab10e] KLU v0.6.0
⌃ [63c18a36] KernelAbstractions v0.9.22
⌃ [ba0b0d4f] Krylov v0.9.6
⌅ [929cbde3] LLVM v8.0.0
⌃ [b964fa9f] LaTeXStrings v1.3.1
[2ee39098] LabelledArrays v1.16.0
⌅ [984bce1d] LambertW v0.4.6
⌃ [23fbe1c1] Latexify v0.16.4
[10f19ff3] LayoutPointers v0.1.17
⌃ [5078a376] LazyArrays v2.1.9
[2d8b4e74] LevyArea v1.0.0
⌃ [d3d80556] LineSearches v7.2.0
⌅ [7ed4a6bd] LinearSolve v2.30.2
⌃ [2ab3a3ac] LogExpFunctions v0.3.28
⌃ [e6f89c97] LoggingExtras v1.0.3
⌃ [bdcacae8] LoopVectorization v0.12.171
⌅ [6f1432cf] LoweredCodeUtils v2.4.8
[d8e11817] MLStyle v0.4.17
⌃ [1914dd2f] MacroTools v0.5.13
[d125e4d3] ManualMemory v0.1.8
⌃ [bb5d69b7] MaybeInplace v0.1.3
[739be429] MbedTLS v1.1.9
[442fdcdd] Measures v0.3.2
[e1d29d7a] Missings v1.2.0
⌅ [dbe65cb8] MistyClosures v1.0.1
⌃ [961ee093] ModelingToolkit v9.26.0
[46d2c3a1] MuladdMacro v0.2.4
⌃ [102ac46a] MultivariatePolynomials v0.5.6
⌃ [ffc61752] Mustache v1.0.19
⌃ [d8a4904e] MutableArithmetics v1.4.5
⌃ [d41bc354] NLSolversBase v7.8.3
[2774e3e8] NLsolve v4.5.1
⌃ [872c559c] NNlib v0.9.21
⌃ [77ba4419] NaNMath v1.0.2
⌅ [8913a72c] NonlinearSolve v3.13.1
⌃ [d8793406] ObjectFile v0.4.1
⌃ [6fe1bfb0] OffsetArrays v1.14.1
[4d8831e6] OpenSSL v1.4.3
⌃ [429524aa] Optim v1.9.4
⌅ [3bd65402] Optimisers v0.3.3
⌃ [bac558e1] OrderedCollections v1.6.3
⌃ [1dea7af3] OrdinaryDiffEq v6.86.0
⌃ [90014a1f] PDMats v0.11.31
[65ce6f38] PackageExtensionCompat v1.0.2
⌃ [65888b18] ParameterizedFunctions v5.17.0
[d96e819e] Parameters v0.12.3
[69de0a69] Parsers v2.8.1
[b98c9c47] Pipe v1.3.0
⌃ [ccf2f8ad] PlotThemes v3.2.0
⌃ [995b91a9] PlotUtils v1.4.1
⌃ [91a5bcdd] Plots v1.40.5
[e409e4f3] PoissonRandom v0.4.4
⌃ [f517fe37] Polyester v0.7.15
[1d0040c9] PolyesterWeave v0.2.2
[2dfb63ee] PooledArrays v1.4.3
[85a6dd25] PositiveFactorizations v0.2.4
⌃ [d236fae5] PreallocationTools v0.4.22
⌅ [aea7be01] PrecompileTools v1.2.1
[21216c6a] Preferences v1.4.3
⌃ [08abe8d2] PrettyTables v2.3.2
[92933f4c] ProgressMeter v1.10.2
⌃ [43287f4e] PtrArrays v1.2.0
⌃ [1fd47b50] QuadGK v2.9.4
[74087812] Random123 v1.7.0
⌃ [e6cf234a] RandomNumbers v1.5.3
[c1ae055f] RealDot v0.1.0
[3cdcf5f2] RecipesBase v1.3.4
[01d81517] RecipesPipeline v0.6.12
⌃ [731186ca] RecursiveArrayTools v3.26.0
[f2c3362d] RecursiveFactorization v0.2.23
[189a3867] Reexport v1.2.2
[05181044] RelocatableFolders v1.0.1
⌃ [ae029012] Requires v1.3.0
[ae5879a3] ResettableStacks v1.1.1
[37e2e3b7] ReverseDiff v1.15.3
⌅ [79098fc4] Rmath v0.7.1
[7e49a35a] RuntimeGeneratedFunctions v0.5.13
[94e857df] SIMDTypes v0.1.0
[476501e8] SLEEFPirates v0.6.43
⌃ [0bca4576] SciMLBase v2.44.0
[31c91b34] SciMLBenchmarks v0.1.3
⌃ [c0aeaf25] SciMLOperators v0.3.8
⌃ [1ed8b502] SciMLSensitivity v7.64.0
⌃ [53ae85a6] SciMLStructures v1.4.1
[6c6a2e73] Scratch v1.2.1
⌃ [91c51154] SentinelArrays v1.4.5
⌃ [efcf1570] Setfield v1.1.1
[992d4aef] Showoff v1.0.3
⌃ [777ac1f9] SimpleBufferStream v1.1.0
⌅ [727e6d20] SimpleNonlinearSolve v1.11.0
[699a6c99] SimpleTraits v0.9.4
[ce78b400] SimpleUnPack v1.1.0
[b85f4697] SoftGlobalScope v1.1.0
[a2af1166] SortingAlgorithms v1.2.1
⌃ [47a9eef4] SparseDiffTools v2.19.0
[dc90abb0] SparseInverseSubset v0.1.2
⌅ [0a514795] SparseMatrixColorings v0.3.5
[e56a9233] Sparspak v0.3.9
⌃ [276daf66] SpecialFunctions v2.4.0
⌃ [aedffcd0] Static v1.1.1
⌃ [0d7ed370] StaticArrayInterface v1.5.1
⌃ [90137ffa] StaticArrays v1.9.7
[1e83bf80] StaticArraysCore v1.4.3
[82ae8749] StatsAPI v1.7.0
⌃ [2913bbd2] StatsBase v0.34.3
⌃ [4c63d2b9] StatsFuns v1.3.1
⌃ [789caeaf] StochasticDiffEq v6.66.0
[7792a7ef] StrideArraysCore v0.5.7
[69024149] StringEncodings v0.3.7
⌅ [892a3eda] StringManipulation v0.3.4
⌅ [09ab397b] StructArrays v0.6.18
⌃ [53d494c1] StructIO v0.3.0
⌃ [2efcf032] SymbolicIndexingInterface v0.3.26
⌃ [19f23fe9] SymbolicLimits v0.2.1
⌅ [d1185830] SymbolicUtils v2.1.2
⌅ [0c5d862f] Symbolics v5.34.0
[9ce81f87] TableMetadataTools v0.1.0
[3783bdb8] TableTraits v1.0.1
[bd369af6] Tables v1.12.0
⌃ [07d77754] Tapir v0.2.26
[62fd8b95] TensorCore v0.1.1
⌅ [8ea1fca8] TermInterface v0.4.1
[8290d209] ThreadingUtilities v0.5.2
⌃ [a759f4b9] TimerOutputs v0.5.24
[0796e94c] Tokenize v0.5.29
⌃ [9f7883ad] Tracker v0.2.34
⌃ [3bb67fe8] TranscodingStreams v0.11.1
[d5829a12] TriangularSolve v0.2.1
⌃ [410a4b4d] Tricks v0.1.8
[781d530d] TruncatedStacktraces v1.4.0
[5c2747f8] URIs v1.5.1
[3a884ed6] UnPack v1.0.2
[1cfade01] UnicodeFun v0.4.1
⌃ [1986cc42] Unitful v1.21.0
[45397f5d] UnitfulLatexify v1.6.4
[a7c27f48] Unityper v0.1.6
⌅ [013be700] UnsafeAtomics v0.2.1
⌅ [d80eeb9a] UnsafeAtomicsLLVM v0.1.5
[41fe7b60] Unzip v0.2.0
⌃ [3d5dd08c] VectorizationBase v0.21.70
[81def892] VersionParsing v1.3.0
[19fa3120] VertexSafeGraphs v0.2.0
[44d3d7a6] Weave v0.10.12
⌃ [ddb6d928] YAML v0.4.11
⌃ [c2297ded] ZMQ v1.2.6
⌅ [e88e6eb3] Zygote v0.6.70
⌃ [700de1a5] ZygoteRules v0.2.5
⌃ [6e34b625] Bzip2_jll v1.0.8+1
⌃ [83423d85] Cairo_jll v1.18.0+2
⌅ [7cc45869] Enzyme_jll v0.0.137+0
⌃ [2702e6a9] EpollShim_jll v0.0.20230411+0
⌃ [2e619515] Expat_jll v2.6.2+0
⌅ [b22a6f82] FFMPEG_jll v4.4.4+1
⌃ [a3f928ae] Fontconfig_jll v2.13.96+0
⌃ [d7e528f0] FreeType2_jll v2.13.2+0
⌃ [559328eb] FriBidi_jll v1.0.14+0
⌃ [0656b61e] GLFW_jll v3.4.0+0
⌅ [d2c73de3] GR_jll v0.73.7+0
[78b55507] Gettext_jll v0.21.0+0
⌃ [f8c6e375] Git_jll v2.44.0+2
⌃ [7746bdde] Glib_jll v2.80.2+0
⌃ [3b182d85] Graphite2_jll v1.3.14+0
⌅ [2e76f6c2] HarfBuzz_jll v2.8.1+1
⌅ [1d5cc7b8] IntelOpenMP_jll v2024.2.0+0
⌃ [aacddb02] JpegTurbo_jll v3.0.3+0
[c1c5ebd0] LAME_jll v3.100.2+0
⌅ [88015f11] LERC_jll v3.0.0+1
⌅ [dad2f222] LLVMExtra_jll v0.0.30+0
⌃ [1d63c593] LLVMOpenMP_jll v15.0.7+0
⌃ [dd4b983a] LZO_jll v2.10.2+0
⌅ [e9f186c6] Libffi_jll v3.2.2+1
⌃ [d4300ac3] Libgcrypt_jll v1.8.11+0
⌃ [7e76a0d4] Libglvnd_jll v1.6.0+0
⌃ [7add5ba3] Libgpg_error_jll v1.49.0+0
⌃ [94ce4f54] Libiconv_jll v1.17.0+0
⌃ [4b2f31a3] Libmount_jll v2.40.1+0
⌅ [89763e89] Libtiff_jll v4.5.1+1
⌃ [38a345b3] Libuuid_jll v2.40.1+0
⌃ [856f044c] MKL_jll v2024.2.0+0
[e7412a2a] Ogg_jll v1.3.5+1
⌃ [458c3c95] OpenSSL_jll v3.0.14+0
⌃ [efe28fd5] OpenSpecFun_jll v0.5.5+0
⌃ [91d4177d] Opus_jll v1.3.2+0
⌅ [30392449] Pixman_jll v0.43.4+0
⌅ [c0090381] Qt6Base_jll v6.7.1+1
⌅ [629bc702] Qt6Declarative_jll v6.7.1+2
⌅ [ce943373] Qt6ShaderTools_jll v6.7.1+1
⌃ [e99dba38] Qt6Wayland_jll v6.7.1+1
⌅ [f50d1b31] Rmath_jll v0.4.2+0
[a44049a8] Vulkan_Loader_jll v1.3.243+0
⌃ [a2964d1f] Wayland_jll v1.21.0+1
⌃ [2381bf8a] Wayland_protocols_jll v1.31.0+0
⌃ [02c8fc9c] XML2_jll v2.13.1+0
⌃ [aed1982a] XSLT_jll v1.1.41+0
⌃ [ffd25f8a] XZ_jll v5.4.6+0
[f67eecfb] Xorg_libICE_jll v1.1.1+0
[c834827a] Xorg_libSM_jll v1.2.4+0
⌃ [4f6342f7] Xorg_libX11_jll v1.8.6+0
⌃ [0c0b7dd1] Xorg_libXau_jll v1.0.11+0
⌃ [935fb764] Xorg_libXcursor_jll v1.2.0+4
⌃ [a3789734] Xorg_libXdmcp_jll v1.1.4+0
⌃ [1082639a] Xorg_libXext_jll v1.3.6+0
⌃ [d091e8ba] Xorg_libXfixes_jll v5.0.3+4
⌃ [a51aa0fd] Xorg_libXi_jll v1.7.10+4
⌃ [d1454406] Xorg_libXinerama_jll v1.1.4+4
⌃ [ec84b674] Xorg_libXrandr_jll v1.5.2+4
⌃ [ea2f1a96] Xorg_libXrender_jll v0.9.11+0
⌃ [14d82f49] Xorg_libpthread_stubs_jll v0.1.1+0
⌃ [c7cfdc94] Xorg_libxcb_jll v1.17.0+0
⌃ [cc61e674] Xorg_libxkbfile_jll v1.1.2+0
[e920d4aa] Xorg_xcb_util_cursor_jll v0.1.4+0
[12413925] Xorg_xcb_util_image_jll v0.4.0+1
[2def613f] Xorg_xcb_util_jll v0.4.0+1
[975044d2] Xorg_xcb_util_keysyms_jll v0.4.0+1
[0d47668e] Xorg_xcb_util_renderutil_jll v0.3.9+1
[c22f9ab0] Xorg_xcb_util_wm_jll v0.4.1+1
⌃ [35661453] Xorg_xkbcomp_jll v1.4.6+0
[33bec58e] Xorg_xkeyboard_config_jll v2.39.0+0
⌃ [c5fb5394] Xorg_xtrans_jll v1.5.0+0
⌃ [8f1865be] ZeroMQ_jll v4.3.5+0
⌃ [3161d3a3] Zstd_jll v1.5.6+0
[35ca27e7] eudev_jll v3.2.9+0
⌅ [214eeab7] fzf_jll v0.43.0+0
⌃ [1a1c6b14] gperf_jll v3.1.1+0
⌃ [a4ae2306] libaom_jll v3.9.0+0
⌃ [0ac62f75] libass_jll v0.15.1+0
[2db6ffa8] libevdev_jll v1.11.0+0
⌃ [f638f0a6] libfdk_aac_jll v2.0.2+0
[36db933b] libinput_jll v1.18.0+0
⌃ [b53b4c65] libpng_jll v1.6.43+1
⌃ [a9144af2] libsodium_jll v1.0.20+0
⌃ [f27f6e37] libvorbis_jll v1.3.7+1
[009596ad] mtdev_jll v1.1.6+0
⌃ [1317d2d5] oneTBB_jll v2021.12.0+0
⌅ [1270edf5] x264_jll v2021.5.5+0
⌅ [dfaa095f] x265_jll v3.5.0+0
⌃ [d8fb68d0] xkbcommon_jll v1.4.1+1
[0dad84c5] ArgTools v1.1.1
[56f22d72] Artifacts
[2a0f44e3] Base64
[ade2ca70] Dates
[8ba89e20] Distributed
[f43a241f] Downloads v1.6.0
[7b1f6079] FileWatching
[9fa8497b] Future
[b77e0a4c] InteractiveUtils
[4af54fe1] LazyArtifacts
[b27032c2] LibCURL v0.6.4
[76f85450] LibGit2
[8f399da3] Libdl
[37e2e46d] LinearAlgebra
[56ddb016] Logging
[d6f4376e] Markdown
[a63ad114] Mmap
[ca575930] NetworkOptions v1.2.0
[44cfe95a] Pkg v1.10.0
[de0858da] Printf
[9abbd945] Profile
[3fa0cd96] REPL
[9a3f8284] Random
[ea8e919c] SHA v0.7.0
[9e88b42a] Serialization
[1a1011a3] SharedArrays
[6462fe0b] Sockets
[2f01184e] SparseArrays v1.10.0
[10745b16] Statistics v1.10.0
[4607b0f0] SuiteSparse
[fa267f1f] TOML v1.0.3
[a4e569a6] Tar v1.10.0
[8dfed614] Test
[cf7118a7] UUIDs
[4ec0a83e] Unicode
[e66e0078] CompilerSupportLibraries_jll v1.1.1+0
[deac9b47] LibCURL_jll v8.4.0+0
[e37daf67] LibGit2_jll v1.6.4+0
[29816b5a] LibSSH2_jll v1.11.0+1
[c8ffd9c3] MbedTLS_jll v2.28.2+1
[14a3606d] MozillaCACerts_jll v2023.1.10
[4536629a] OpenBLAS_jll v0.3.23+4
[05823500] OpenLibm_jll v0.8.1+2
[efcefdf7] PCRE2_jll v10.42.0+1
[bea87d4a] SuiteSparse_jll v7.2.1+1
[83775a58] Zlib_jll v1.2.13+1
[8e850b90] libblastrampoline_jll v5.8.0+1
[8e850ede] nghttp2_jll v1.52.0+1
[3f19e933] p7zip_jll v17.4.0+2
Info Packages marked with ⌃ and ⌅ have new versions available. Those with ⌃ may be upgradable, but those with ⌅ are restricted by compatibility constraints from upgrading. To see why use `status --outdated -m`
Warning The project dependencies or compat requirements have changed since the manifest was last resolved. It is recommended to `Pkg.resolve()` or consider `Pkg.update()` if necessary.