Brusselator Work-Precision Diagrams
using OrdinaryDiffEq, DiffEqDevTools, Sundials, ParameterizedFunctions, Plots,
ODEInterfaceDiffEq, LSODA, SparseArrays, LinearSolve,
LinearAlgebra, IncompleteLU, AlgebraicMultigrid, Symbolics, ModelingToolkit,
RecursiveFactorization
gr()
const N = 8
xyd_brusselator = range(0,stop=1,length=N)
brusselator_f(x, y, t) = (((x-0.3)^2 + (y-0.6)^2) <= 0.1^2) * (t >= 1.1) * 5.
limit(a, N) = a == N+1 ? 1 : a == 0 ? N : a
function brusselator_2d_loop(du, u, p, t)
A, B, alpha, dx = p
alpha = alpha/dx^2
@inbounds for I in CartesianIndices((N, N))
i, j = Tuple(I)
x, y = xyd_brusselator[I[1]], xyd_brusselator[I[2]]
ip1, im1, jp1, jm1 = limit(i+1, N), limit(i-1, N), limit(j+1, N), limit(j-1, N)
du[i,j,1] = alpha*(u[im1,j,1] + u[ip1,j,1] + u[i,jp1,1] + u[i,jm1,1] - 4u[i,j,1]) +
B + u[i,j,1]^2*u[i,j,2] - (A + 1)*u[i,j,1] + brusselator_f(x, y, t)
du[i,j,2] = alpha*(u[im1,j,2] + u[ip1,j,2] + u[i,jp1,2] + u[i,jm1,2] - 4u[i,j,2]) +
A*u[i,j,1] - u[i,j,1]^2*u[i,j,2]
end
end
p = (3.4, 1., 10., step(xyd_brusselator))
input = rand(N,N,2)
output = similar(input)
sparsity_pattern = Symbolics.jacobian_sparsity(brusselator_2d_loop,output,input,p,0.0)
jac_sparsity = Float64.(sparse(sparsity_pattern))
f = ODEFunction{true, SciMLBase.FullSpecialize}(brusselator_2d_loop;jac_prototype=jac_sparsity)
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
u
end
u0 = init_brusselator_2d(xyd_brusselator)
prob = ODEProblem(f,u0,(0.,11.5),p);
prob_mtk = ODEProblem(complete(modelingtoolkitize(prob)),[],(0.0,11.5),jac=true,sparse=true);
Also comparing with MethodOfLines.jl:
using MethodOfLines, DomainSets
@parameters x y t
@variables u(..) v(..)
Dt = Differential(t)
Dx = Differential(x)
Dy = Differential(y)
Dxx = Differential(x)^2
Dyy = Differential(y)^2
∇²(u) = Dxx(u) + Dyy(u)
brusselator_f(x, y, t) = (((x-0.3)^2 + (y-0.6)^2) <= 0.1^2) * (t >= 1.1) * 5.
x_min = y_min = t_min = 0.0
x_max = y_max = 1.0
t_max = 11.5
α = 10.
u0_mol(x,y,t) = 22(y*(1-y))^(3/2)
v0_mol(x,y,t) = 27(x*(1-x))^(3/2)
eq = [Dt(u(x,y,t)) ~ 1. + v(x,y,t)*u(x,y,t)^2 - 4.4*u(x,y,t) + α*∇²(u(x,y,t)) + brusselator_f(x, y, t),
Dt(v(x,y,t)) ~ 3.4*u(x,y,t) - v(x,y,t)*u(x,y,t)^2 + α*∇²(v(x,y,t))]
domains = [x ∈ Interval(x_min, x_max),
y ∈ Interval(y_min, y_max),
t ∈ Interval(t_min, t_max)]
bcs = [u(x,y,0) ~ u0_mol(x,y,0),
u(0,y,t) ~ u(1,y,t),
u(x,0,t) ~ u(x,1,t),
v(x,y,0) ~ v0_mol(x,y,0),
v(0,y,t) ~ v(1,y,t),
v(x,0,t) ~ v(x,1,t)]
@named pdesys = PDESystem(eq,bcs,domains,[x,y,t],[u(x,y,t),v(x,y,t)])
# Method of lines discretization
dx = 1/N
dy = 1/N
order = 2
discretization = MOLFiniteDifference([x=>dx, y=>dy], t; approx_order = order, jac = true, sparse = true, wrap = Val(false))
# Convert the PDE system into an ODE problem
prob_mol = discretize(pdesys,discretization)
ODEProblem with uType Vector{Float64} and tType Float64. In-place: true
timespan: (0.0, 11.5)
u0: 128-element Vector{Float64}:
0.7957923865311465
0.7957923865311465
0.7957923865311465
0.7957923865311465
0.7957923865311465
0.7957923865311465
0.7957923865311465
0.7957923865311465
1.7861773953054048
1.7861773953054048
⋮
0.0
0.9766542925609525
2.1921268033293604
3.063590342214851
3.3750000000000004
3.063590342214851
2.1921268033293604
0.9766542925609525
0.0
using Base.Experimental: Const, @aliasscope
macro vp(expr)
nodes = (Symbol("llvm.loop.vectorize.predicate.enable"), 1)
if expr.head != :for
error("Syntax error: loopinfo needs a for loop")
end
push!(expr.args[2].args, Expr(:loopinfo, nodes))
return esc(expr)
end
struct Brusselator2DLoop <: Function
N::Int
s::Float64
end
function (b::Brusselator2DLoop)(du, unc, p, t)
N = b.N
s = b.s
A, B, alpha, dx = p
alpha = alpha/abs2(dx)
u = Base.Experimental.Const(unc)
Base.Experimental.@aliasscope begin
@inbounds @fastmath begin
b = ((abs2(-0.3) + abs2(-0.6)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du1 = alpha*(u[N,1,1] + u[2,1,1] + u[1,2,1] + u[1,N,1] - 4u[1,1,1]) +
B + abs2(u[1,1,1])*u[1,1,2] - (A + 1)*u[1,1,1] + b
du2 = alpha*(u[N,1,2] + u[2,1,2] + u[1,2,2] + u[1,N,2] - 4u[1,1,2]) +
A*u[1,1,1] - abs2(u[1,1,1])*u[1,1,2]
du[1,1,1] = du1
du[1,1,2] = du2
@vp for i = 2:N-1
x = (i-1)*s
ip1 = i+1
im1 = i-1
b = ((abs2(x-0.3) + abs2(-0.6)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du1 = alpha*(u[im1,1,1] + u[ip1,1,1] + u[i,2,1] + u[i,N,1] - 4u[i,1,1]) +
B + abs2(u[i,1,1])*u[i,1,2] - (A + 1)*u[i,1,1] + b
du2 = alpha*(u[im1,1,2] + u[ip1,1,2] + u[i,2,2] + u[i,N,2] - 4u[i,1,2]) +
A*u[i,1,1] - abs2(u[i,1,1])*u[i,1,2]
du[i,1,1] = du1
du[i,1,2] = du2
end
b = ((abs2(0.7) + abs2(-0.6)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du1 = alpha*(u[N-1,1,1] + u[1,1,1] + u[N,2,1] + u[N,N,1] - 4u[N,1,1]) +
B + abs2(u[N,1,1])*u[N,1,2] - (A + 1)*u[N,1,1] + b
du2 = alpha*(u[N-1,1,2] + u[1,1,2] + u[N,2,2] + u[N,N,2] - 4u[N,1,2]) +
A*u[N,1,1] - abs2(u[N,1,1])*u[N,1,2]
du[N,1,1] = du1
du[N,1,2] = du2
for j = 2:N-1
y = (j-1)*s
jp1 = j+1
jm1 = j-1
b0 = ((abs2(-0.3) + abs2(y-0.6)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du[1,j,1] = alpha*(u[N,j,1] + u[2,j,1] + u[1,jp1,1] + u[1,jm1,1] - 4u[1,j,1]) +
B + abs2(u[1,j,1])*u[1,j,2] - (A + 1)*u[1,j,1] + b0
du[1,j,2] = alpha*(u[N,j,2] + u[2,j,2] + u[1,jp1,2] + u[1,jm1,2] - 4u[1,j,2]) +
A*u[1,j,1] - abs2(u[1,j,1])*u[1,j,2]
@vp for i = 2:N-1
x = (i-1)*s
b = ((abs2(x-0.3) + abs2(y-0.6)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du1 = alpha*(u[i-1,j,1] + u[i+1,j,1] + u[i,jp1,1] + u[i,jm1,1] - 4u[i,j,1]) +
B + abs2(u[i,j,1])*u[i,j,2] - (A + 1)*u[i,j,1] + b
du2 = alpha*(u[i-1,j,2] + u[i+1,j,2] + u[i,jp1,2] + u[i,jm1,2] - 4u[i,j,2]) +
A*u[i,j,1] - abs2(u[i,j,1])*u[i,j,2]
du[i,j,1] = du1
du[i,j,2] = du2
end
bN = ((abs2(0.7) + abs2(y-0.6)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du[N,j,1] = alpha*(u[N-1,j,1] + u[1,j,1] + u[N,jp1,1] + u[N,jm1,1] - 4u[N,j,1]) +
B + abs2(u[N,j,1])*u[N,j,2] - (A + 1)*u[N,j,1] + bN
du[N,j,2] = alpha*(u[N-1,j,2] + u[1,j,2] + u[N,jp1,2] + u[N,jm1,2] - 4u[N,j,2]) +
A*u[N,j,1] - abs2(u[N,j,1])*u[N,j,2]
end
b = ((abs2(-0.3) + abs2(0.4)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du1 = alpha*(u[N,N,1] + u[2,N,1] + u[1,1,1] + u[1,N-1,1] - 4u[1,N,1]) +
B + abs2(u[1,N,1])*u[1,N,2] - (A + 1)*u[1,N,1] + b
du2 = alpha*(u[N,N,2] + u[2,N,2] + u[1,1,2] + u[1,N-1,2] - 4u[1,N,2]) +
A*u[1,N,1] - abs2(u[1,N,1])*u[1,N,2]
du[1,N,1] = du1
du[1,N,2] = du2
@vp for i = 2:N-1
x = (i-1)*s
ip1 = i+1
im1 = i-1
b = ((abs2(x-0.3) + abs2(0.4)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du1 = alpha*(u[im1,N,1] + u[ip1,N,1] + u[i,1,1] + u[i,N-1,1] - 4u[i,N,1]) +
B + abs2(u[i,N,1])*u[i,N,2] - (A + 1)*u[i,N,1] + b
du2 = alpha*(u[im1,N,2] + u[ip1,N,2] + u[i,1,2] + u[i,N-1,2] - 4u[i,N,2]) +
A*u[i,N,1] - abs2(u[i,N,1])*u[i,N,2]
du[i,N,1] = du1
du[i,N,2] = du2
end
b = ((abs2(0.7) + abs2(0.4)) <= abs2(0.1)) * (t >= 1.1) * 5.0
du1 = alpha*(u[N-1,N,1] + u[1,N,1] + u[N,1,1] + u[N,N-1,1] - 4u[N,N,1]) +
B + abs2(u[N,N,1])*u[N,N,2] - (A + 1)*u[N,N,1] + b
du2 = alpha*(u[N-1,N,2] + u[1,N,2] + u[N,1,2] + u[N,N-1,2] - 4u[N,N,2]) +
A*u[N,N,1] - abs2(u[N,N,1])*u[N,N,2]
du[N,N,1] = du1
du[N,N,2] = du2
end
end
end
function fast_bruss(N)
xyd_brusselator = range(0,stop=1,length=N)
brusselator_2d_loop = Brusselator2DLoop(N,Float64(step(xyd_brusselator)))
p = (3.4, 1., 10., step(xyd_brusselator))
input = rand(N,N,2)
output = similar(input)
sparsity_pattern = Symbolics.jacobian_sparsity(brusselator_2d_loop,output,input,p,0.0)
jac_sparsity = Float64.(sparse(sparsity_pattern))
f = ODEFunction(brusselator_2d_loop;jac_prototype=jac_sparsity)
u0 = zeros(N, N, 2)
@inbounds for I in CartesianIndices((N, N))
x = xyd_brusselator[I[1]]
y = xyd_brusselator[I[2]]
u0[I,1] = 22*(y*(1-y))^(3/2)
u0[I,2] = 27*(x*(1-x))^(3/2)
end
return ODEProblem(f,u0,(0.,11.5),p)
end
fastprob = fast_bruss(N)
ODEProblem with uType Array{Float64, 3} and tType Float64. In-place: true
timespan: (0.0, 11.5)
u0: 8×8×2 Array{Float64, 3}:
[:, :, 1] =
0.0 0.942661 2.02828 2.66625 2.66625 2.02828 0.942661 0.0
0.0 0.942661 2.02828 2.66625 2.66625 2.02828 0.942661 0.0
0.0 0.942661 2.02828 2.66625 2.66625 2.02828 0.942661 0.0
0.0 0.942661 2.02828 2.66625 2.66625 2.02828 0.942661 0.0
0.0 0.942661 2.02828 2.66625 2.66625 2.02828 0.942661 0.0
0.0 0.942661 2.02828 2.66625 2.66625 2.02828 0.942661 0.0
0.0 0.942661 2.02828 2.66625 2.66625 2.02828 0.942661 0.0
0.0 0.942661 2.02828 2.66625 2.66625 2.02828 0.942661 0.0
[:, :, 2] =
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
1.1569 1.1569 1.1569 1.1569 1.1569 1.1569 1.1569 1.1569
2.48926 2.48926 2.48926 2.48926 2.48926 2.48926 2.48926 2.48926
3.27221 3.27221 3.27221 3.27221 3.27221 3.27221 3.27221 3.27221
3.27221 3.27221 3.27221 3.27221 3.27221 3.27221 3.27221 3.27221
2.48926 2.48926 2.48926 2.48926 2.48926 2.48926 2.48926 2.48926
1.1569 1.1569 1.1569 1.1569 1.1569 1.1569 1.1569 1.1569
0.0 0.0 0.0 0.0 0.0 0.0 0.0 0.0
sol = solve(prob,CVODE_BDF(),abstol=1/10^14,reltol=1/10^14)
sol2 = solve(prob_mtk,CVODE_BDF(linear_solver = :KLU),abstol=1/10^14,reltol=1/10^14)
sol3 = solve(prob_mol,CVODE_BDF(linear_solver = :KLU),abstol=1/10^14,reltol=1/10^14,wrap=Val(false))
retcode: Success
Interpolation: 3rd order Hermite
t: 7088-element Vector{Float64}:
0.0
1.4238870862469186e-15
1.4240294749555432e-11
1.5662900337424728e-10
2.9901771199893913e-10
5.160801983264086e-10
8.725711849204022e-10
1.4781897844618872e-9
2.519675048784877e-9
4.3190349346343956e-9
⋮
11.495686108181351
11.496274749829004
11.496863391476657
11.49745203312431
11.498040674771962
11.498629316419615
11.499217958067268
11.49980659971492
11.5
u: 7088-element Vector{Vector{Float64}}:
[0.7957923865311465, 0.7957923865311465, 0.7957923865311465, 0.79579238653
11465, 0.7957923865311465, 0.7957923865311465, 0.7957923865311465, 0.795792
3865311465, 1.7861773953054048, 1.7861773953054048 … 0.9766542925609525,
0.0, 0.9766542925609525, 2.1921268033293604, 3.063590342214851, 3.375000000
0000004, 3.063590342214851, 2.1921268033293604, 0.9766542925609525, 0.0]
[0.7957923865313211, 0.7957923865313222, 0.795792386531323, 0.795792386531
3233, 0.795792386531323, 0.7957923865313222, 0.7957923865313211, 0.79579238
65313202, 1.786177395305144, 1.7861773953051494 … 0.9766542925611731, 1.7
838787595843077e-12, 0.9766542925611701, 2.192126803329047, 3.0635903422143
41, 3.374999999999433, 3.063590342214341, 2.192126803329047, 0.976654292561
1701, 1.7800261566757715e-12]
[0.7957923882778083, 0.7957923882887697, 0.7957923882966287, 0.79579238829
94371, 0.7957923882966287, 0.7957923882887697, 0.7957923882778083, 0.795792
3882690007, 1.7861773926974722, 1.7861773927526943 … 0.9766542947672135,
1.7840571189693174e-8, 0.9766542947374913, 2.1921268001941354, 3.0635903371
106386, 3.374999994323757, 3.0635903371106386, 2.1921268001941354, 0.976654
2947374913, 1.7802041309561312e-8]
[0.7957924057426853, 0.7957924058632492, 0.7957924059496904, 0.79579240598
05794, 0.7957924059496904, 0.7957924058632492, 0.7957924057426853, 0.795792
4056458101, 1.786177366620756, 1.7861773672281442 … 0.9766543168276247, 1
.9622842530506597e-7, 0.9766543165007099, 2.1921267688450223, 3.06359028607
36164, 3.3749999375669977, 3.0635902860736164, 2.1921267688450223, 0.976654
3165007099, 1.958046352534055e-7]
[0.795792423207569, 0.7957924234377355, 0.7957924236027589, 0.795792423661
7286, 0.7957924236027589, 0.7957924234377355, 0.795792423207569, 0.79579242
30226263, 1.7861773405440422, 1.7861773417035964 … 0.9766543388880443, 3.
746162540951588e-7, 0.976654338263937, 2.192126737495912, 3.063590235036595
6, 3.3749998808102393, 3.0635902350365956, 2.192126737495912, 0.97665433826
3937, 3.738072040103176e-7]
[0.7957924498316827, 0.7957924502289312, 0.7957924505137483, 0.79579245061
55252, 0.7957924505137483, 0.7957924502289312, 0.7957924498316827, 0.795792
4495124866, 1.7861773007917632, 1.7861773027930594 … 0.9766543725177618,
6.465570447995342e-7, 0.9766543714406037, 2.1921266897062024, 3.06359015723
39215, 3.37499979428819, 3.0635901572339215, 2.1921266897062024, 0.97665437
14406037, 6.451606907088272e-7]
[0.7957924935576395, 0.795792494229294, 0.7957924947108534, 0.795792494882
9345, 0.7957924947108534, 0.795792494229294, 0.7957924935576395, 0.79579249
30179533, 1.7861772355049113, 1.7861772388886359 … 0.9766544277493292, 1.
0931768703959916e-6, 0.9766544259281074, 2.1921266112191384, 3.063590029455
275, 3.3749996521893415, 3.063590029455275, 2.1921266112191384, 0.976654425
9281074, 1.0908159626472516e-6]
[0.7957925678408412, 0.7957925689786661, 0.7957925697944583, 0.79579257008
59744, 0.7957925697944583, 0.7957925689786661, 0.7957925678408412, 0.795792
5669265793, 1.7861771245934985, 1.786177130325737 … 0.9766545215786651, 1
.8519088939775023e-6, 0.9766545184934051, 2.192126477882781, 3.063589812380
733, 3.3749994107871655, 3.063589812380733, 2.192126477882781, 0.9766545184
934051, 1.8479093730188316e-6]
[0.7957926955863502, 0.7957926975258519, 0.7957926989164267, 0.79579269941
33361, 0.7957926989164267, 0.7957926975258519, 0.7957926955863502, 0.795792
6940279286, 1.7861769338587068, 1.786176943629694 … 0.976654682937829, 3.
1567028571118565e-6, 0.9766546776788025, 2.1921262485837283, 3.063589439076
6394, 3.374998995646712, 3.0635894390766394, 2.1921262485837283, 0.97665467
76788025, 3.149885411590933e-6]
[0.7957929162914903, 0.795792919616041, 0.7957929219996626, 0.795792922851
4283, 0.7957929219996626, 0.795792919616041, 0.7957929162914903, 0.79579291
36201628, 1.7861766043291272, 1.786176621077798 … 0.9766549617170381, 5.4
10973999659674e-6, 0.9766549527024425, 2.19212585242728, 3.0635887941244206
, 3.374998278414355, 3.0635887941244206, 2.19212585242728, 0.97665495270244
25, 5.3992880839237515e-6]
⋮
[4.073542147517943, 4.073678256826177, 4.073678256826177, 4.07354214751794
25, 4.073369276100977, 4.073256699954195, 4.073256699954196, 4.073369276100
977, 4.073684991585881, 4.073870247312897 … 1.1851911298275786, 1.1851854
979836303, 1.185186253669876, 1.1851814266888845, 1.1851814266888845, 1.185
186253669876, 1.1851926200032414, 1.185196887133102, 1.185196887133102, 1.1
851926200032414]
[4.075236586856852, 4.075372679859942, 4.075372679859942, 4.07523658685685
2, 4.075063737468067, 4.07495117632641, 4.074951176326411, 4.07506373746806
7, 4.075379416824593, 4.075564653507394 … 1.1817785913918286, 1.181772981
8020342, 1.1817737345428994, 1.1817689266686777, 1.1817689266686775, 1.1817
737345428994, 1.1817800756316128, 1.1817843258181333, 1.1817843258181333, 1
.1817800756316128]
[4.076902965031779, 4.077039041822214, 4.077039041822214, 4.07690296503177
85, 4.076730137546133, 4.076617591324378, 4.076617591324378, 4.076730137546
133, 4.077045780979023, 4.077230998725429 … 1.1783931250209605, 1.1783875
375594715, 1.1783882873714766, 1.1783834984958645, 1.1783834984958643, 1.17
83882873714766, 1.1783946033585904, 1.178398836697857, 1.178398836697857, 1
.1783946033585904]
[4.0785414099893, 4.078677470659699, 4.078677470659699, 4.078541409989299,
4.078368604281576, 4.078256072894378, 4.078256072894378, 4.078368604281576
, 4.078684211995861, 4.078869410913847 … 1.1750346192486743, 1.1750290537
894597, 1.1750298006891513, 1.1750250307038326, 1.1750250307038328, 1.17502
98006891513, 1.1750360917179217, 1.1750403083061565, 1.1750403083061562, 1.
1750360917179214]
[4.0801520507468645, 4.080288095389964, 4.080288095389964, 4.0801520507468
64, 4.079979266691691, 4.079866750053601, 4.079866750053601, 4.079979266691
691, 4.080294838892661, 4.080480019090341 … 1.1717029614621586, 1.1716974
178790238, 1.1716981618829718, 1.1716934106794916, 1.1716934106794918, 1.17
16981618829718, 1.1717044280968363, 1.1717086280303821, 1.171708628030382,
1.1717044280968363]
[4.081735017358575, 4.0818710460672145, 4.0818710460672145, 4.081735017358
575, 4.081562254830442, 4.081449752855921, 4.081449752855921, 4.08156225483
0442, 4.0818777917236195, 4.082062953309229 … 1.168398037935694, 1.168392
5161023, 1.168393257227095, 1.1683885246968762, 1.1683885246968762, 1.16839
32572270952, 1.1683994987696518, 1.168403682144958, 1.168403682144958, 1.16
83994987696518]
[4.08329044088123, 4.083426453748336, 4.083426453748336, 4.08329044088123,
4.083117699754509, 4.083005212357935, 4.083005212357935, 4.083117699754509
, 4.0834332015456125, 4.083618344627493 … 1.1651197338640091, 1.165114233
653892, 1.1651149719161427, 1.165110257950501, 1.165110257950501, 1.1651149
719161427, 1.165121188931128, 1.1651253558447368, 1.165125355844737, 1.1651
211889311281]
[4.084818453340643, 4.0849544504592155, 4.084954450459216, 4.0848184533406
42, 4.084645733489605, 4.08453326058529, 4.08453326058529, 4.08464573348960
6, 4.084961200384522, 4.085146325071104 … 1.1618679333953803, 1.161862454
6819691, 1.1618631900983, 1.1618584945884618, 1.1618584945884618, 1.1618631
900982999, 1.161869382729568, 1.1618735332780994, 1.1618735332780996, 1.161
869382729568]
[4.085314527232731, 4.085450519197376, 4.0854505191973765, 4.0853145272327
3, 4.085141814344439, 4.085029346182895, 4.085029346182895, 4.0851418143444
4, 4.0854572698190506, 4.0856423884854705 … 1.1608053075612759, 1.1607998
358830436, 1.160800570367968, 1.160795880898096, 1.160795880898096, 1.16080
0570367968, 1.1608067550192933, 1.1608109002121196, 1.1608109002121196, 1.1
608067550192935]
test_sol = [sol,sol2,sol,sol3]
probs = [prob,prob_mtk,fastprob,prob_mol];
plot(sol, idxs = 1)
plot(sol, idxs = 10)
Setup Preconditioners
OrdinaryDiffEq
function incompletelu(W,du,u,p,t,newW,Plprev,Prprev,solverdata)
if newW === nothing || newW
Pl = ilu(convert(AbstractMatrix,W), τ = 50.0)
else
Pl = Plprev
end
Pl,nothing
end
function algebraicmultigrid(W,du,u,p,t,newW,Plprev,Prprev,solverdata)
if newW === nothing || newW
Pl = aspreconditioner(ruge_stuben(convert(AbstractMatrix,W)))
else
Pl = Plprev
end
Pl,nothing
end
algebraicmultigrid (generic function with 1 method)
Sundials
const jaccache = prob_mtk.f.jac(prob.u0,prob.p,0.0)
const W = I - 1.0*jaccache
prectmp = ilu(W, τ = 50.0)
const preccache = Ref(prectmp)
function psetupilu(p, t, u, du, jok, jcurPtr, gamma)
if !jok
prob_mtk.f.jac(jaccache,u,p,t)
jcurPtr[] = true
# W = I - gamma*J
@. W = -gamma*jaccache
idxs = diagind(W)
@. @view(W[idxs]) = @view(W[idxs]) + 1
# Build preconditioner on W
preccache[] = ilu(W, τ = 5.0)
end
end
function precilu(z,r,p,t,y,fy,gamma,delta,lr)
ldiv!(z,preccache[],r)
end
prectmp2 = aspreconditioner(ruge_stuben(W, presmoother = AlgebraicMultigrid.Jacobi(rand(size(W,1))), postsmoother = AlgebraicMultigrid.Jacobi(rand(size(W,1)))))
const preccache2 = Ref(prectmp2)
function psetupamg(p, t, u, du, jok, jcurPtr, gamma)
if !jok
prob_mtk.f.jac(jaccache,u,p,t)
jcurPtr[] = true
# W = I - gamma*J
@. W = -gamma*jaccache
idxs = diagind(W)
@. @view(W[idxs]) = @view(W[idxs]) + 1
# Build preconditioner on W
preccache2[] = aspreconditioner(ruge_stuben(W, presmoother = AlgebraicMultigrid.Jacobi(rand(size(W,1))), postsmoother = AlgebraicMultigrid.Jacobi(rand(size(W,1)))))
end
end
function precamg(z,r,p,t,y,fy,gamma,delta,lr)
ldiv!(z,preccache2[],r)
end
precamg (generic function with 1 method)
Compare Problem Implementations
abstols = 1.0 ./ 10.0 .^ (5:8)
reltols = 1.0 ./ 10.0 .^ (1:4);
setups = [
Dict(:alg => KenCarp47(linsolve=KLUFactorization())),
Dict(:alg => KenCarp47(linsolve=KLUFactorization()), :prob_choice => 2),
Dict(:alg => KenCarp47(linsolve=KLUFactorization()), :prob_choice => 3),
Dict(:alg => KenCarp47(linsolve=KLUFactorization()), :prob_choice => 4),
Dict(:alg => KenCarp47(linsolve=KrylovJL_GMRES())),
Dict(:alg => KenCarp47(linsolve=KrylovJL_GMRES()), :prob_choice => 2),
Dict(:alg => KenCarp47(linsolve=KrylovJL_GMRES()), :prob_choice => 3),
Dict(:alg => KenCarp47(linsolve=KrylovJL_GMRES()), :prob_choice => 4),]
names = ["KenCarp47 KLU","KenCarp47 KLU MTK","KenCarp47 KLU FastBruss", "KenCarp47 KLU MOL",
"KenCarp47 GMRES", "KenCarp47 GMRES MTK", "KenCarp47 GMRES FastBruss", "KenCarp47 GMRES MOL"];
wp = WorkPrecisionSet(probs,abstols,reltols,setups;names = names,
save_everystep=false,appxsol=test_sol,maxiters=Int(1e5),numruns=10,wrap=Val(false))
plot(wp)
High Tolerances
This is the speed when you just want the answer.
abstols = 1.0 ./ 10.0 .^ (5:8)
reltols = 1.0 ./ 10.0 .^ (1:4);
setups = [
Dict(:alg=>CVODE_BDF(linear_solver = :KLU), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver = :GMRES)),
Dict(:alg=>CVODE_BDF(linear_solver = :GMRES), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precilu,psetup=psetupilu,prec_side=1)),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precamg,psetup=psetupamg,prec_side=1)),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precilu,psetup=psetupilu,prec_side=1), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precamg,psetup=psetupamg,prec_side=1), :prob_choice => 2),
]
names = ["CVODE MTK KLU","CVODE GMRES","CVODE MTK GMRES", "CVODE iLU GMRES", "CVODE AMG GMRES", "CVODE iLU MTK GMRES", "CVODE AMG MTK GMRES"];
wp = WorkPrecisionSet(probs,abstols,reltols,setups;names=names,
save_everystep=false,appxsol=test_sol,maxiters=Int(1e5),numruns=10)
plot(wp)
setups = [
Dict(:alg=>KenCarp47(linsolve=KLUFactorization())),
Dict(:alg=>KenCarp47(linsolve=KLUFactorization()), :prob_choice => 2),
Dict(:alg=>KenCarp47(linsolve=UMFPACKFactorization())),
Dict(:alg=>KenCarp47(linsolve=UMFPACKFactorization()), :prob_choice => 2),
Dict(:alg=>KenCarp47(linsolve=KrylovJL_GMRES())),
Dict(:alg=>KenCarp47(linsolve=KrylovJL_GMRES()), :prob_choice => 2),
Dict(:alg=>KenCarp47(linsolve=KrylovJL_GMRES(),precs=incompletelu,concrete_jac=true)),
Dict(:alg=>KenCarp47(linsolve=KrylovJL_GMRES(),precs=incompletelu,concrete_jac=true), :prob_choice => 2),
Dict(:alg=>KenCarp47(linsolve=KrylovJL_GMRES(),precs=algebraicmultigrid,concrete_jac=true)),
Dict(:alg=>KenCarp47(linsolve=KrylovJL_GMRES(),precs=algebraicmultigrid,concrete_jac=true), :prob_choice => 2),
]
names = ["KenCarp47 KLU","KenCarp47 KLU MTK","KenCarp47 UMFPACK", "KenCarp47 UMFPACK MTK", "KenCarp47 GMRES",
"KenCarp47 GMRES MTK", "KenCarp47 iLU GMRES", "KenCarp47 iLU GMRES MTK", "KenCarp47 AMG GMRES",
"KenCarp47 AMG GMRES MTK"];
wp = WorkPrecisionSet(probs,abstols,reltols,setups;names = names,
save_everystep=false,appxsol=test_sol,maxiters=Int(1e5),numruns=10)
plot(wp)
setups = [
Dict(:alg=>TRBDF2()),
Dict(:alg=>KenCarp4()),
Dict(:alg=>KenCarp47()),
# Dict(:alg=>QNDF()), # bad
Dict(:alg=>FBDF()),
]
wp = WorkPrecisionSet(probs,abstols,reltols,setups;
save_everystep=false,appxsol=test_sol,maxiters=Int(1e5),numruns=10)
plot(wp)
setups = [
Dict(:alg=>KenCarp47(linsolve=KLUFactorization()), :prob_choice => 2),
Dict(:alg=>KenCarp47(linsolve=KrylovJL_GMRES()), :prob_choice => 2),
Dict(:alg=>FBDF(linsolve=KLUFactorization()), :prob_choice => 2),
Dict(:alg=>FBDF(linsolve=KrylovJL_GMRES()), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver = :KLU), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precilu,psetup=psetupilu,prec_side=1), :prob_choice => 2),
]
names = ["KenCarp47 KLU MTK", "KenCarp47 GMRES MTK",
"FBDF KLU MTK", "FBDF GMRES MTK",
"CVODE MTK KLU", "CVODE iLU MTK GMRES"
];
wp = WorkPrecisionSet(probs,abstols,reltols,setups;names = names,
save_everystep=false,appxsol=test_sol,maxiters=Int(1e5),numruns=10)
plot(wp)
Low Tolerances
This is the speed at lower tolerances, measuring what's good when accuracy is needed.
abstols = 1.0 ./ 10.0 .^ (7:12)
reltols = 1.0 ./ 10.0 .^ (4:9)
setups = [
Dict(:alg=>CVODE_BDF(linear_solver = :KLU), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver = :GMRES)),
Dict(:alg=>CVODE_BDF(linear_solver = :GMRES), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precilu,psetup=psetupilu,prec_side=1)),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precamg,psetup=psetupamg,prec_side=1)),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precilu,psetup=psetupilu,prec_side=1), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precamg,psetup=psetupamg,prec_side=1), :prob_choice => 2),
]
names = ["CVODE MTK KLU","CVODE GMRES","CVODE MTK GMRES", "CVODE iLU GMRES", "CVODE AMG GMRES", "CVODE iLU MTK GMRES", "CVODE AMG MTK GMRES"];
wp = WorkPrecisionSet(probs,abstols,reltols,setups;names = names,
save_everystep=false,appxsol=test_sol,maxiters=Int(1e5),numruns=10)
plot(wp)
setups = [
Dict(:alg=>KenCarp47(linsolve=KLUFactorization()), :prob_choice => 2),
Dict(:alg=>KenCarp47(linsolve=KrylovJL_GMRES()), :prob_choice => 2),
Dict(:alg=>FBDF(linsolve=KLUFactorization()), :prob_choice => 2),
Dict(:alg=>FBDF(linsolve=KrylovJL_GMRES()), :prob_choice => 2),
Dict(:alg=>Rodas5P(linsolve=KrylovJL_GMRES()), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver = :KLU), :prob_choice => 2),
Dict(:alg=>CVODE_BDF(linear_solver=:GMRES,prec=precilu,psetup=psetupilu,prec_side=1), :prob_choice => 2),
]
names = ["KenCarp47 KLU MTK", "KenCarp47 GMRES MTK",
"FBDF KLU MTK", "FBDF GMRES MTK",
"Rodas5P GMRES MTK",
"CVODE MTK KLU", "CVODE iLU MTK GMRES"
];
wp = WorkPrecisionSet(probs,abstols,reltols,setups;names = names,
save_everystep=false,appxsol=test_sol,maxiters=Int(1e5),numruns=10)
plot(wp)
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/StiffODE","Bruss.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/StiffODE/Project.toml`
⌃ [2169fc97] AlgebraicMultigrid v0.6.0
⌃ [6e4b80f9] BenchmarkTools v1.5.0
⌃ [f3b72e0c] DiffEqDevTools v2.45.0
⌃ [5b8099bc] DomainSets v0.7.14
[5a33fad7] GeometricIntegratorsDiffEq v0.2.5
[40713840] IncompleteLU v0.2.1
[7f56f5a3] LSODA v0.7.5
⌅ [7ed4a6bd] LinearSolve v2.34.0
⌃ [94925ecb] MethodOfLines v0.11.3
⌃ [961ee093] ModelingToolkit v9.32.0
⌃ [09606e27] ODEInterfaceDiffEq v3.13.3
⌃ [1dea7af3] OrdinaryDiffEq v6.89.0
⌃ [65888b18] ParameterizedFunctions v5.17.0
⌃ [91a5bcdd] Plots v1.40.8
[132c30aa] ProfileSVG v0.2.2
[f2c3362d] RecursiveFactorization v0.2.23
[31c91b34] SciMLBenchmarks v0.1.3
⌃ [90137ffa] StaticArrays v1.9.7
⌃ [c3572dad] Sundials v4.25.0
⌅ [0c5d862f] Symbolics v5.36.0
⌃ [a759f4b9] TimerOutputs v0.5.24
[37e2e46d] LinearAlgebra
[2f01184e] SparseArrays v1.10.0
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/StiffODE/Manifest.toml`
⌃ [47edcb42] ADTypes v1.7.1
[a4c015fc] ANSIColoredPrinters v0.0.1
[621f4979] AbstractFFTs v1.5.0
[1520ce14] AbstractTrees v0.4.5
⌃ [7d9f7c33] Accessors v0.1.37
⌃ [79e6a3ab] Adapt v4.0.4
⌃ [2169fc97] AlgebraicMultigrid v0.6.0
[66dad0bd] AliasTables v1.1.3
⌅ [4c88cf16] Aqua v0.5.6
[ec485272] ArnoldiMethod v0.4.0
⌃ [4fba245c] ArrayInterface v7.16.0
⌃ [4c555306] ArrayLayouts v1.10.3
[13072b0f] AxisAlgorithms v1.1.0
⌃ [aae01518] BandedMatrices v1.7.4
⌃ [6e4b80f9] BenchmarkTools v1.5.0
[e2ed5e7c] Bijections v0.1.9
[d1d4a3ce] BitFlags v0.1.9
[62783981] BitTwiddlingConvenienceFunctions v0.1.6
⌃ [8e7c35d0] BlockArrays v1.1.0
[fa961155] CEnum v0.5.0
[2a0fbf3d] CPUSummary v0.2.6
[00ebfdb7] CSTParser v3.4.3
⌃ [d360d2e6] ChainRulesCore v1.24.0
[fb6a15b2] CloseOpenIntervals v0.1.13
⌃ [944b1d66] CodecZlib v0.7.6
⌃ [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.1
[f70d9fcc] CommonWorldInvalidations v1.0.0
⌃ [a09551c4] CompactBasisFunctions v0.2.12
[34da2185] Compat v4.16.0
[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
⌅ [7ae1f121] ContinuumArrays v0.18.4
[d38c429a] Contour v0.6.3
[adafc99b] CpuId v0.3.1
[a8cc5b0e] Crayons v4.1.1
⌅ [717857b8] DSP v0.7.9
[9a962f9c] DataAPI v1.16.0
⌃ [864edb3b] DataStructures v0.18.20
[e2d170a0] DataValueInterfaces v1.0.0
⌃ [55939f99] DecFP v1.3.2
[8bb1440f] DelimitedFiles v1.9.1
⌃ [2b5f629d] DiffEqBase v6.155.0
⌅ [459566f4] DiffEqCallbacks v3.9.1
⌃ [f3b72e0c] DiffEqDevTools v2.45.0
⌃ [77a26b50] DiffEqNoiseProcess v5.23.0
[163ba53b] DiffResults v1.1.0
[b552c78f] DiffRules v1.15.1
⌅ [a0c0ee7d] DifferentiationInterface v0.5.17
⌃ [b4f34e82] Distances v0.10.11
⌃ [31c24e10] Distributions v0.25.111
[ffbed154] DocStringExtensions v0.9.3
⌅ [e30172f5] Documenter v0.27.25
⌃ [5b8099bc] DomainSets v0.7.14
⌅ [7c1d4256] DynamicPolynomials v0.5.7
⌅ [06fc5a27] DynamicQuantities v0.13.2
[4e289a0a] EnumX v1.0.4
⌅ [f151be2c] EnzymeCore v0.7.8
⌃ [460bff9d] ExceptionUnwrapping v0.1.10
⌃ [d4d017d3] ExponentialUtilities v1.26.1
[e2ba6199] ExprTools v0.1.10
⌅ [6b7a57c9] Expronicon v0.8.5
⌃ [c87230d0] FFMPEG v0.4.1
⌃ [7a1cc6ca] FFTW v1.8.0
[7034ab61] FastBroadcast v0.3.5
[9aa1b823] FastClosures v0.3.2
[442a2c76] FastGaussQuadrature v1.0.2
[29a986be] FastLapackInterface v2.0.4
⌅ [057dd010] FastTransforms v0.16.4
⌃ [5789e2e9] FileIO v1.16.3
[1a297f60] FillArrays v1.13.0
[64ca27bc] FindFirstFunctions v1.4.1
⌃ [6a86dc24] FiniteDiff v2.24.0
[53c48c17] FixedPointNumbers v0.8.5
⌃ [08572546] FlameGraphs v1.0.0
[1fa38f19] Format v1.3.7
⌃ [f6369f11] ForwardDiff v0.10.36
[069b7b12] FunctionWrappers v1.1.3
[77dc65aa] FunctionWrappersWrappers v0.1.3
⌅ [d9f16b24] Functors v0.4.12
⌅ [46192b85] GPUArraysCore v0.1.6
⌃ [28b8d3ca] GR v0.73.7
[a8297547] GenericFFT v0.1.6
⌃ [14197337] GenericLinearAlgebra v0.3.12
[c145ed77] GenericSchur v0.5.4
⌅ [9a0b12b7] GeometricBase v0.2.6
⌅ [c85262ba] GeometricEquations v0.2.1
⌅ [dcce2d33] GeometricIntegrators v0.9.2
[5a33fad7] GeometricIntegratorsDiffEq v0.2.5
[d7ba0133] Git v1.3.1
[c27321d9] Glob v1.3.1
⌃ [86223c79] Graphs v1.11.2
[42e2da0e] Grisu v1.0.2
⌅ [f67ccb44] HDF5 v0.16.16
⌃ [cd3eb016] HTTP v1.10.8
[eafb193a] Highlights v0.5.3
[3e5b6fbb] HostCPUFeatures v0.1.17
⌃ [34004b35] HypergeometricFunctions v0.3.24
⌃ [7073ff75] IJulia v1.25.0
[b5f81e59] IOCapture v0.2.5
[615f187c] IfElse v0.1.1
[40713840] IncompleteLU v0.2.1
[9b13fd28] IndirectArrays v1.0.0
⌅ [4858937d] InfiniteArrays v0.14.3
⌃ [e1ba4f0e] Infinities v0.1.8
[d25df0c9] Inflate v0.1.5
[a98d9a8b] Interpolations v0.15.1
[8197267c] IntervalSets v0.7.10
⌃ [3587e190] InverseFunctions v0.1.16
⌃ [92d709cd] IrrationalConstants v0.2.2
[c8e1da08] IterTools v1.10.0
[82899510] IteratorInterfaceExtensions v1.0.0
⌃ [1019f520] JLFzf v0.1.8
⌃ [692b3bcd] JLLWrappers v1.6.0
[682c06a0] JSON v0.21.4
⌃ [98e50ef6] JuliaFormatter v1.0.60
⌃ [ccbc3e58] JumpProcesses v9.13.7
[ef3ab10e] KLU v0.6.0
⌃ [ba0b0d4f] Krylov v0.9.6
[7f56f5a3] LSODA v0.7.5
⌃ [b964fa9f] LaTeXStrings v1.3.1
[2ee39098] LabelledArrays v1.16.0
⌅ [984bce1d] LambertW v0.4.6
⌃ [23fbe1c1] Latexify v0.16.5
[10f19ff3] LayoutPointers v0.1.17
⌃ [5078a376] LazyArrays v2.2.1
[1d6d02ad] LeftChildRightSiblingTrees v0.2.0
[d3d80556] LineSearches v7.3.0
⌅ [7ed4a6bd] LinearSolve v2.34.0
⌃ [2ab3a3ac] LogExpFunctions v0.3.28
⌃ [e6f89c97] LoggingExtras v1.0.3
[bdcacae8] LoopVectorization v0.12.171
[d8e11817] MLStyle v0.4.17
[3da0fdf6] MPIPreferences v0.1.11
⌃ [1914dd2f] MacroTools v0.5.13
[d125e4d3] ManualMemory v0.1.8
[bb5d69b7] MaybeInplace v0.1.4
[739be429] MbedTLS v1.1.9
[442fdcdd] Measures v0.3.2
⌃ [94925ecb] MethodOfLines v0.11.3
[e1d29d7a] Missings v1.2.0
⌃ [961ee093] ModelingToolkit v9.32.0
[46d2c3a1] MuladdMacro v0.2.4
⌃ [102ac46a] MultivariatePolynomials v0.5.6
[ffc61752] Mustache v1.0.20
⌃ [d8a4904e] MutableArithmetics v1.4.6
⌃ [d41bc354] NLSolversBase v7.8.3
[2774e3e8] NLsolve v4.5.1
⌃ [77ba4419] NaNMath v1.0.2
⌅ [8913a72c] NonlinearSolve v3.14.0
[54ca160b] ODEInterface v0.5.0
⌃ [09606e27] ODEInterfaceDiffEq v3.13.3
⌃ [6fe1bfb0] OffsetArrays v1.14.1
[4d8831e6] OpenSSL v1.4.3
⌃ [429524aa] Optim v1.9.4
⌃ [bac558e1] OrderedCollections v1.6.3
⌃ [1dea7af3] OrdinaryDiffEq v6.89.0
⌃ [89bda076] OrdinaryDiffEqAdamsBashforthMoulton v1.1.0
⌃ [6ad6398a] OrdinaryDiffEqBDF v1.1.2
⌃ [bbf590c4] OrdinaryDiffEqCore v1.6.0
⌃ [50262376] OrdinaryDiffEqDefault v1.1.0
⌃ [4302a76b] OrdinaryDiffEqDifferentiation v1.1.0
[9286f039] OrdinaryDiffEqExplicitRK v1.1.0
⌃ [e0540318] OrdinaryDiffEqExponentialRK v1.1.0
⌃ [becaefa8] OrdinaryDiffEqExtrapolation v1.1.0
⌃ [5960d6e9] OrdinaryDiffEqFIRK v1.1.1
[101fe9f7] OrdinaryDiffEqFeagin v1.1.0
[d3585ca7] OrdinaryDiffEqFunctionMap v1.1.1
[d28bc4f8] OrdinaryDiffEqHighOrderRK v1.1.0
⌃ [9f002381] OrdinaryDiffEqIMEXMultistep v1.1.0
[521117fe] OrdinaryDiffEqLinear v1.1.0
[1344f307] OrdinaryDiffEqLowOrderRK v1.2.0
⌃ [b0944070] OrdinaryDiffEqLowStorageRK v1.2.1
⌃ [127b3ac7] OrdinaryDiffEqNonlinearSolve v1.2.1
[c9986a66] OrdinaryDiffEqNordsieck v1.1.0
⌃ [5dd0a6cf] OrdinaryDiffEqPDIRK v1.1.0
[5b33eab2] OrdinaryDiffEqPRK v1.1.0
[04162be5] OrdinaryDiffEqQPRK v1.1.0
[af6ede74] OrdinaryDiffEqRKN v1.1.0
⌃ [43230ef6] OrdinaryDiffEqRosenbrock v1.2.0
⌃ [2d112036] OrdinaryDiffEqSDIRK v1.1.0
⌃ [669c94d9] OrdinaryDiffEqSSPRK v1.2.0
⌃ [e3e12d00] OrdinaryDiffEqStabilizedIRK v1.1.0
[358294b1] OrdinaryDiffEqStabilizedRK v1.1.0
⌃ [fa646aed] OrdinaryDiffEqSymplecticRK v1.1.0
[b1df2697] OrdinaryDiffEqTsit5 v1.1.0
[79d7bb75] OrdinaryDiffEqVerner v1.1.1
⌃ [a7812802] PDEBase v0.1.14
⌃ [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.8
[e409e4f3] PoissonRandom v0.4.4
[f517fe37] Polyester v0.7.16
[1d0040c9] PolyesterWeave v0.2.2
⌅ [f27b6e38] Polynomials v3.2.13
[85a6dd25] PositiveFactorizations v0.2.4
⌃ [d236fae5] PreallocationTools v0.4.24
[aea7be01] PrecompileTools v1.2.1
[21216c6a] Preferences v1.4.3
⌃ [08abe8d2] PrettyTables v2.3.2
[132c30aa] ProfileSVG v0.2.2
[92933f4c] ProgressMeter v1.10.2
⌃ [43287f4e] PtrArrays v1.2.1
⌃ [1fd47b50] QuadGK v2.11.0
[a08977f5] QuadratureRules v0.1.6
⌅ [c4ea9172] QuasiArrays v0.11.7
[74087812] Random123 v1.7.0
[e6cf234a] RandomNumbers v1.6.0
[c84ed2f1] Ratios v0.4.5
[3cdcf5f2] RecipesBase v1.3.4
[01d81517] RecipesPipeline v0.6.12
⌃ [731186ca] RecursiveArrayTools v3.27.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
⌅ [79098fc4] Rmath v0.7.1
[47965b36] RootedTrees v2.23.1
⌅ [fb486d5c] RungeKutta v0.4.6
[7e49a35a] RuntimeGeneratedFunctions v0.5.13
[94e857df] SIMDTypes v0.1.0
[476501e8] SLEEFPirates v0.6.43
⌃ [0bca4576] SciMLBase v2.53.0
[31c91b34] SciMLBenchmarks v0.1.3
⌃ [c0aeaf25] SciMLOperators v0.3.10
⌃ [53ae85a6] SciMLStructures v1.5.0
[6c6a2e73] Scratch v1.2.1
⌃ [efcf1570] Setfield v1.1.1
[992d4aef] Showoff v1.0.3
⌃ [777ac1f9] SimpleBufferStream v1.1.0
⌅ [727e6d20] SimpleNonlinearSolve v1.12.0
⌅ [36b790f5] SimpleSolvers v0.2.4
[699a6c99] SimpleTraits v0.9.4
[ce78b400] SimpleUnPack v1.1.0
[b85f4697] SoftGlobalScope v1.1.0
[a2af1166] SortingAlgorithms v1.2.1
⌃ [47a9eef4] SparseDiffTools v2.20.0
⌃ [0a514795] SparseMatrixColorings v0.4.0
[e56a9233] Sparspak v0.3.9
⌃ [276daf66] SpecialFunctions v2.4.0
⌃ [aedffcd0] Static v1.1.1
[0d7ed370] StaticArrayInterface v1.8.0
⌃ [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
[7792a7ef] StrideArraysCore v0.5.7
[69024149] StringEncodings v0.3.7
⌅ [892a3eda] StringManipulation v0.3.4
⌅ [09ab397b] StructArrays v0.6.18
⌃ [c3572dad] Sundials v4.25.0
⌃ [2efcf032] SymbolicIndexingInterface v0.3.30
⌃ [19f23fe9] SymbolicLimits v0.2.1
⌅ [d1185830] SymbolicUtils v2.1.3
⌅ [0c5d862f] Symbolics v5.36.0
[3783bdb8] TableTraits v1.0.1
[bd369af6] Tables v1.12.0
[62fd8b95] TensorCore v0.1.1
⌅ [8ea1fca8] TermInterface v0.4.1
[8290d209] ThreadingUtilities v0.5.2
⌃ [a759f4b9] TimerOutputs v0.5.24
⌃ [c751599d] ToeplitzMatrices v0.8.4
[0796e94c] Tokenize v0.5.29
⌃ [3bb67fe8] TranscodingStreams v0.11.2
[d5829a12] TriangularSolve v0.2.1
⌃ [410a4b4d] Tricks v0.1.9
[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
[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
[efce3f68] WoodburyMatrices v1.0.0
⌃ [ddb6d928] YAML v0.4.12
⌃ [c2297ded] ZMQ v1.3.0
⌃ [6e34b625] Bzip2_jll v1.0.8+1
⌃ [83423d85] Cairo_jll v1.18.0+2
[ee1fde0b] Dbus_jll v1.14.10+0
[47200ebd] DecFP_jll v2.0.3+1
⌃ [2702e6a9] EpollShim_jll v0.0.20230411+0
⌃ [2e619515] Expat_jll v2.6.2+0
⌅ [b22a6f82] FFMPEG_jll v4.4.4+1
⌃ [f5851436] FFTW_jll v3.3.10+0
[34b6f7d7] FastTransforms_jll v0.6.3+0
⌃ [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+1
⌅ [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
⌃ [0234f1f7] HDF5_jll v1.14.2+1
⌃ [2e76f6c2] HarfBuzz_jll v8.3.1+0
⌃ [e33a78d0] Hwloc_jll v2.11.1+0
⌅ [1d5cc7b8] IntelOpenMP_jll v2024.2.1+0
⌃ [aacddb02] JpegTurbo_jll v3.0.3+0
[c1c5ebd0] LAME_jll v3.100.2+0
⌅ [88015f11] LERC_jll v3.0.0+1
⌃ [1d63c593] LLVMOpenMP_jll v17.0.6+0
[aae0fff6] LSODA_jll v0.1.2+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
⌃ [7cb0a576] MPICH_jll v4.2.2+0
⌃ [f1f71cc9] MPItrampoline_jll v5.4.0+0
⌃ [9237b28f] MicrosoftMPI_jll v10.1.4+2
[c771fb93] ODEInterface_jll v0.0.1+0
[e7412a2a] Ogg_jll v1.3.5+1
⌃ [fe0851c0] OpenMPI_jll v5.0.5+0
⌃ [458c3c95] OpenSSL_jll v3.0.15+0
⌃ [efe28fd5] OpenSpecFun_jll v0.5.5+0
[91d4177d] Opus_jll v1.3.3+0
⌃ [36c8627f] Pango_jll v1.54.1+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.3+0
⌅ [fb77eaff] Sundials_jll v5.2.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.3+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.53.0+0
⌃ [1a1c6b14] gperf_jll v3.1.1+0
⌃ [477f73a3] libaec_jll v1.1.2+0
⌃ [a4ae2306] libaom_jll v3.9.0+0
[0ac62f75] libass_jll v0.15.2+0
[1183f4f0] libdecor_jll v0.2.2+0
[2db6ffa8] libevdev_jll v1.11.0+0
[f638f0a6] libfdk_aac_jll v2.0.3+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+2
[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
[781609d7] GMP_jll v6.2.1+6
[deac9b47] LibCURL_jll v8.4.0+0
[e37daf67] LibGit2_jll v1.6.4+0
[29816b5a] LibSSH2_jll v1.11.0+1
[3a97d323] MPFR_jll v4.2.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.