Brusselator Work-Precision Diagrams

using OrdinaryDiffEq, DiffEqDevTools, Sundials, ParameterizedFunctions, Plots,
      ODEInterfaceDiffEq, LSODA, SparseArrays, LinearSolve,
      LinearAlgebra, IncompleteLU, AlgebraicMultigrid, Symbolics, ModelingToolkit
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(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: 7100-element Vector{Float64}:
  0.0
  1.4238870862469186e-15
  1.4240294749555432e-11
  1.5662900337424728e-10
  2.9901771199893913e-10
  5.160801983204045e-10
  8.725711850742434e-10
  1.4781897823229784e-9
  2.5196750517210653e-9
  4.3190349441260325e-9
  ⋮
 11.495766698882717
 11.496356788238053
 11.49694687759339
 11.497536966948726
 11.498127056304062
 11.498717145659398
 11.499307235014735
 11.499897324370071
 11.5
u: 7100-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.465570447920121e-7, 0.9766543714406037, 2.1921266897062024, 3.06359015723
39215, 3.37499979428819, 3.0635901572339215, 2.1921266897062024, 0.97665437
14406037, 6.451606907013214e-7]
 [0.7957924935576395, 0.7957924942292941, 0.7957924947108534, 0.79579249488
29346, 0.7957924947108534, 0.7957924942292941, 0.7957924935576395, 0.795792
4930179534, 1.786177235504911, 1.7861772388886359  …  0.9766544277493292, 1
.0931768705887272e-6, 0.9766544259281074, 2.1921266112191384, 3.06359002945
5275, 3.3749996521893415, 3.063590029455275, 2.1921266112191384, 0.97665442
59281074, 1.0908159628395709e-6]
 [0.795792567840841, 0.795792568978666, 0.795792569794458, 0.79579257008597
43, 0.795792569794458, 0.795792568978666, 0.795792567840841, 0.795792566926
5791, 1.7861771245934988, 1.7861771303257374  …  0.9766545215786647, 1.8519
088912978322e-6, 0.9766545184934048, 2.1921264778827814, 3.0635898123807337
, 3.3749994107871664, 3.0635898123807337, 2.1921264778827814, 0.97665451849
34048, 1.8479093703449485e-6]
 [0.7957926955863507, 0.7957926975258524, 0.795792698916427, 0.795792699413
3366, 0.795792698916427, 0.7957926975258524, 0.7957926955863507, 0.79579269
40279291, 1.7861769338587061, 1.7861769436296935  …  0.9766546829378294, 3.
156702860790371e-6, 0.976654677678803, 2.192126248583728, 3.063589439076638
5, 3.3749989956467106, 3.0635894390766385, 2.192126248583728, 0.97665467767
8803, 3.149885415261503e-6]
 [0.7957929162914915, 0.7957929196160424, 0.7957929219996638, 0.79579292285
14298, 0.7957929219996638, 0.7957929196160424, 0.7957929162914915, 0.795792
913620164, 1.7861766043291254, 1.7861766210777965  …  0.9766549617170396, 5
.410974011550959e-6, 0.976654952702444, 2.192125852427278, 3.06358879412441
75, 3.374998278414351, 3.0635887941244175, 2.192125852427278, 0.97665495270
2444, 5.399288095789355e-6]
 ⋮
 [4.073775795169845, 4.0739119022402885, 4.0739119022402885, 4.073775795169
845, 4.07360292677613, 4.073490352688726, 4.073490352688725, 4.073602926776
13, 4.073918637302626, 4.074103890415941  …  1.1847223169188434, 1.18471668
81291311, 1.1847174434111514, 1.1847126190524446, 1.1847126190524444, 1.184
7174434111511, 1.1847238062798366, 1.1847280710843193, 1.1847280710843193, 
1.1847238062798369]
 [4.075470508367493, 4.07560659910549, 4.07560659910549, 4.075470508367493,
 4.075297662038831, 4.075185102981661, 4.075185102981661, 4.075297662038831
, 4.075613336376429, 4.075798570413612  …  1.1813051411720268, 1.1812995346
73785, 1.1813002870054663, 1.1812954817855528, 1.1812954817855525, 1.181300
287005466, 1.18130662458721, 1.1813108724199792, 1.1813108724199792, 1.1813
066245872101]
 [4.07713704005818, 4.077273114556696, 4.077273114556696, 4.07713704005818,
 4.076964215668887, 4.07685167155629, 4.07685167155629, 4.0769642156688874,
 4.077279854023378, 4.0774650690928045  …  1.1779151552745553, 1.1779095709
412886, 1.1779103203392371, 1.1779055341494007, 1.1779055341494007, 1.17791
0320339237, 1.1779166327778257, 1.1779208637354608, 1.1779208637354608, 1.1
779166327778257]
 [4.078775519288571, 4.078911577640698, 4.078911577640698, 4.07877551928857
1, 4.0786027167127905, 4.078490187458987, 4.078490187458987, 4.078602716712
7905, 4.0789183192902545, 4.079103515500452  …  1.1745522467712588, 1.17454
66844762906, 1.1745474309571382, 1.1745426636885088, 1.1745426636885088, 1.
1745474309571382, 1.1745537183965602, 1.1745579325757753, 1.174557932575775
3, 1.1745537183965602]
 [4.080386076181861, 4.0805221184808085, 4.0805221184808085, 4.080386076181
8616, 4.080213295293583, 4.08010078081269, 4.08010078081269, 4.080213295293
582, 4.080528862300358, 4.08071403975999  …  1.1712163020539625, 1.17121076
1670455, 1.1712115052508563, 1.171206756794426, 1.171206756794426, 1.171211
5052508563, 1.1712177678352804, 1.1712219653329095, 1.1712219653329095, 1.1
712177678352806]
 [4.08196884190317, 4.082104868242247, 4.082104868242248, 4.081968841903171
, 4.081796082576247, 4.08168358278229, 4.081683582782289, 4.081796082576247
, 4.082111614218899, 4.082296773036751  …  1.1679072063954714, 1.1679016877
96444, 1.1679024284930737, 1.1678976987397138, 1.1678976987397138, 1.167902
4284930737, 1.1679086663668272, 1.1679128472798093, 1.167912847279809, 1.16
79086663668274]
 [4.083523948625205, 4.083659959097808, 4.083659959097809, 4.08352394862520
55, 4.083351210733373, 4.083238725540298, 4.083238725540297, 4.083351210733
373, 4.083666707218661, 4.083851847503623  …  1.1646248439833002, 1.1646193
470416486, 1.1646200848712, 1.164615373711677, 1.164615373711677, 1.1646200
848711998, 1.164626298178746, 1.164630462604111, 1.1646304626041108, 1.1646
262981787463]
 [4.085051529494204, 4.085187524193804, 4.085187524193804, 4.08505152949420
5, 4.084878812911104, 4.08476634223279, 4.08476634223279, 4.084878812911104
, 4.085194274445952, 4.085379396306998  …  1.1613690979531432, 1.1613636225
416584, 1.16136435752084, 1.1613596648458318, 1.1613596648458318, 1.1613643
5752084, 1.1613705464067576, 1.1613746944416117, 1.1613746944416117, 1.1613
705464067576]
 [4.085314527250261, 4.085450519214907, 4.085450519214907, 4.08531452725026
2, 4.085141814361971, 4.085029346200427, 4.085029346200426, 4.0851418143619
71, 4.085457269836582, 4.085642388503001  …  1.1608053075274634, 1.16079983
5849231, 1.1608005703341557, 1.160795880864284, 1.160795880864284, 1.160800
5703341557, 1.1608067549854808, 1.1608109001783067, 1.1608109001783067, 1.1
608067549854808]
test_sol = [sol,sol2,sol,sol3]
probs = [prob,prob_mtk,fastprob,prob_mol];
plot(sol,vars = 1)

plot(sol,vars = 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.9.3
Commit bed2cd540a1 (2023-08-24 14:43 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-14.0.6 (ORCJIT, znver2)
  Threads: 128 on 128 virtual cores
Environment:
  JULIA_CPU_THREADS = 128
  JULIA_DEPOT_PATH = /cache/julia-buildkite-plugin/depots/5b300254-1738-4989-ae0a-f4d2d937f953
  JULIA_IMAGE_THREADS = 1

Package Information:

Status `/cache/build/exclusive-amdci3-0/julialang/scimlbenchmarks-dot-jl/benchmarks/StiffODE/Project.toml`
⌃ [2169fc97] AlgebraicMultigrid v0.5.1
  [6e4b80f9] BenchmarkTools v1.3.2
⌃ [f3b72e0c] DiffEqDevTools v2.35.0
  [5b8099bc] DomainSets v0.6.7
  [5a33fad7] GeometricIntegratorsDiffEq v0.2.5
  [40713840] IncompleteLU v0.2.1
⌃ [7f56f5a3] LSODA v0.7.4
⌃ [7ed4a6bd] LinearSolve v2.4.2
⌃ [94925ecb] MethodOfLines v0.9.5
⌃ [961ee093] ModelingToolkit v8.64.0
  [c030b06c] ODE v2.15.0
⌃ [09606e27] ODEInterfaceDiffEq v3.13.2
⌃ [1dea7af3] OrdinaryDiffEq v6.53.4
⌃ [65888b18] ParameterizedFunctions v5.15.0
⌃ [91a5bcdd] Plots v1.38.17
  [132c30aa] ProfileSVG v0.2.1
  [31c91b34] SciMLBenchmarks v0.1.3
⌃ [90137ffa] StaticArrays v1.6.2
⌃ [c3572dad] Sundials v4.19.3
⌃ [0c5d862f] Symbolics v5.5.1
  [a759f4b9] TimerOutputs v0.5.23
  [37e2e46d] LinearAlgebra
  [2f01184e] SparseArrays
Info Packages marked with ⌃ have new versions available and may be upgradable.
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-amdci3-0/julialang/scimlbenchmarks-dot-jl/benchmarks/StiffODE/Manifest.toml`
⌅ [47edcb42] ADTypes v0.1.6
  [a4c015fc] ANSIColoredPrinters v0.0.1
⌅ [c3fe647b] AbstractAlgebra v0.31.0
  [621f4979] AbstractFFTs v1.5.0
⌅ [1520ce14] AbstractTrees v0.3.4
  [79e6a3ab] Adapt v3.6.2
⌃ [2169fc97] AlgebraicMultigrid v0.5.1
⌅ [4c88cf16] Aqua v0.5.6
  [ec485272] ArnoldiMethod v0.2.0
  [4fba245c] ArrayInterface v7.4.11
  [30b0a656] ArrayInterfaceCore v0.1.29
⌃ [4c555306] ArrayLayouts v1.1.1
  [13072b0f] AxisAlgorithms v1.0.1
⌃ [aae01518] BandedMatrices v0.17.34
  [6e4b80f9] BenchmarkTools v1.3.2
⌃ [e2ed5e7c] Bijections v0.1.4
  [d1d4a3ce] BitFlags v0.1.7
  [62783981] BitTwiddlingConvenienceFunctions v0.1.5
⌃ [8e7c35d0] BlockArrays v0.16.36
  [fa961155] CEnum v0.4.2
⌃ [2a0fbf3d] CPUSummary v0.2.3
  [00ebfdb7] CSTParser v3.3.6
  [49dc2e85] Calculus v0.5.1
  [d360d2e6] ChainRulesCore v1.16.0
  [fb6a15b2] CloseOpenIntervals v0.1.12
  [944b1d66] CodecZlib v0.7.2
⌃ [35d6a980] ColorSchemes v3.22.0
  [3da002f7] ColorTypes v0.11.4
  [c3611d14] ColorVectorSpace v0.10.0
  [5ae59095] Colors v0.12.10
  [861a8166] Combinatorics v1.0.2
  [a80b9123] CommonMark v0.8.12
  [38540f10] CommonSolve v0.2.4
  [bbf7d656] CommonSubexpressions v0.3.0
⌃ [a09551c4] CompactBasisFunctions v0.2.8
⌃ [34da2185] Compat v4.8.0
  [b152e2b5] CompositeTypes v0.1.3
  [f0e56b4a] ConcurrentUtilities v2.2.1
  [8f4d0f93] Conda v1.9.1
⌃ [187b0558] ConstructionBase v1.5.3
⌅ [7ae1f121] ContinuumArrays v0.14.1
  [d38c429a] Contour v0.6.2
  [adafc99b] CpuId v0.3.1
  [a8cc5b0e] Crayons v4.1.1
⌃ [717857b8] DSP v0.7.8
  [9a962f9c] DataAPI v1.15.0
⌃ [864edb3b] DataStructures v0.18.14
  [e2d170a0] DataValueInterfaces v1.0.0
  [55939f99] DecFP v1.3.2
  [8bb1440f] DelimitedFiles v1.9.1
⌃ [2b5f629d] DiffEqBase v6.127.0
⌃ [459566f4] DiffEqCallbacks v2.27.0
⌃ [f3b72e0c] DiffEqDevTools v2.35.0
⌃ [77a26b50] DiffEqNoiseProcess v5.18.0
  [163ba53b] DiffResults v1.1.0
  [b552c78f] DiffRules v1.15.1
  [b4f34e82] Distances v0.10.9
⌃ [31c24e10] Distributions v0.25.98
  [ffbed154] DocStringExtensions v0.9.3
⌅ [e30172f5] Documenter v0.27.25
  [5b8099bc] DomainSets v0.6.7
  [fa6b7ba4] DualNumbers v0.6.8
⌃ [7c1d4256] DynamicPolynomials v0.5.2
  [4e289a0a] EnumX v1.0.4
  [460bff9d] ExceptionUnwrapping v0.1.9
⌃ [d4d017d3] ExponentialUtilities v1.24.0
  [e2ba6199] ExprTools v0.1.10
  [c87230d0] FFMPEG v0.4.1
  [7a1cc6ca] FFTW v1.7.1
  [7034ab61] FastBroadcast v0.2.6
  [9aa1b823] FastClosures v0.3.2
  [442a2c76] FastGaussQuadrature v0.5.1
  [29a986be] FastLapackInterface v2.0.0
⌃ [057dd010] FastTransforms v0.15.6
  [5789e2e9] FileIO v1.16.1
⌃ [1a297f60] FillArrays v1.5.0
  [6a86dc24] FiniteDiff v2.21.1
  [53c48c17] FixedPointNumbers v0.8.4
⌅ [08572546] FlameGraphs v0.2.10
  [59287772] Formatting v0.4.2
⌃ [f6369f11] ForwardDiff v0.10.35
  [069b7b12] FunctionWrappers v1.1.3
  [77dc65aa] FunctionWrappersWrappers v0.1.3
  [46192b85] GPUArraysCore v0.1.5
⌃ [28b8d3ca] GR v0.72.9
  [a8297547] GenericFFT v0.1.4
  [14197337] GenericLinearAlgebra v0.3.11
  [c145ed77] GenericSchur v0.5.3
⌅ [9a0b12b7] GeometricBase v0.2.6
⌅ [c85262ba] GeometricEquations v0.2.1
⌅ [dcce2d33] GeometricIntegrators v0.9.2
  [5a33fad7] GeometricIntegratorsDiffEq v0.2.5
  [d7ba0133] Git v1.3.0
  [c27321d9] Glob v1.3.1
  [86223c79] Graphs v1.8.0
  [42e2da0e] Grisu v1.0.2
⌃ [0b43b601] Groebner v0.4.2
  [d5909c97] GroupsCore v0.4.0
⌅ [f67ccb44] HDF5 v0.16.15
⌃ [cd3eb016] HTTP v1.9.14
  [eafb193a] Highlights v0.5.2
⌃ [3e5b6fbb] HostCPUFeatures v0.1.15
  [34004b35] HypergeometricFunctions v0.3.23
  [7073ff75] IJulia v1.24.2
  [b5f81e59] IOCapture v0.2.3
  [615f187c] IfElse v0.1.1
  [40713840] IncompleteLU v0.2.1
  [9b13fd28] IndirectArrays v1.0.0
⌃ [4858937d] InfiniteArrays v0.13.0
  [e1ba4f0e] Infinities v0.1.7
  [d25df0c9] Inflate v0.1.3
  [18e54dd8] IntegerMathUtils v0.1.2
  [a98d9a8b] Interpolations v0.14.7
  [8197267c] IntervalSets v0.7.7
  [92d709cd] IrrationalConstants v0.2.2
  [c8e1da08] IterTools v1.8.0
  [82899510] IteratorInterfaceExtensions v1.0.0
  [1019f520] JLFzf v0.1.5
⌃ [692b3bcd] JLLWrappers v1.4.1
  [682c06a0] JSON v0.21.4
⌃ [98e50ef6] JuliaFormatter v1.0.34
⌃ [ccbc3e58] JumpProcesses v9.7.2
⌃ [ef3ab10e] KLU v0.4.0
⌃ [ba0b0d4f] Krylov v0.9.2
⌃ [7f56f5a3] LSODA v0.7.4
  [b964fa9f] LaTeXStrings v1.3.0
  [2ee39098] LabelledArrays v1.14.0
  [984bce1d] LambertW v0.4.6
⌅ [23fbe1c1] Latexify v0.15.21
  [10f19ff3] LayoutPointers v0.1.14
  [50d2b5c4] Lazy v0.15.1
⌃ [5078a376] LazyArrays v1.5.2
⌅ [1d6d02ad] LeftChildRightSiblingTrees v0.1.3
  [d3d80556] LineSearches v7.2.0
⌃ [7ed4a6bd] LinearSolve v2.4.2
⌃ [2ab3a3ac] LogExpFunctions v0.3.24
⌃ [e6f89c97] LoggingExtras v1.0.0
  [bdcacae8] LoopVectorization v0.12.165
  [d8e11817] MLStyle v0.4.17
⌃ [1914dd2f] MacroTools v0.5.10
  [d125e4d3] ManualMemory v0.1.8
⌃ [a3b82374] MatrixFactorizations v2.0.1
  [739be429] MbedTLS v1.1.7
  [442fdcdd] Measures v0.3.2
⌃ [94925ecb] MethodOfLines v0.9.5
  [e1d29d7a] Missings v1.1.0
⌃ [961ee093] ModelingToolkit v8.64.0
  [46d2c3a1] MuladdMacro v0.2.4
⌃ [102ac46a] MultivariatePolynomials v0.5.1
  [ffc61752] Mustache v1.0.17
⌃ [d8a4904e] MutableArithmetics v1.3.0
  [d41bc354] NLSolversBase v7.8.3
  [2774e3e8] NLsolve v4.5.1
  [77ba4419] NaNMath v1.0.2
⌃ [8913a72c] NonlinearSolve v1.9.0
  [c030b06c] ODE v2.15.0
  [54ca160b] ODEInterface v0.5.0
⌃ [09606e27] ODEInterfaceDiffEq v3.13.2
  [6fe1bfb0] OffsetArrays v1.12.10
  [4d8831e6] OpenSSL v1.4.1
⌃ [429524aa] Optim v1.7.6
  [bac558e1] OrderedCollections v1.6.2
⌃ [1dea7af3] OrdinaryDiffEq v6.53.4
  [a7812802] PDEBase v0.1.4
⌃ [90014a1f] PDMats v0.11.17
⌃ [65ce6f38] PackageExtensionCompat v1.0.0
⌃ [65888b18] ParameterizedFunctions v5.15.0
  [d96e819e] Parameters v0.12.3
  [69de0a69] Parsers v2.7.2
  [b98c9c47] Pipe v1.3.0
  [ccf2f8ad] PlotThemes v3.1.0
  [995b91a9] PlotUtils v1.3.5
⌃ [91a5bcdd] Plots v1.38.17
  [e409e4f3] PoissonRandom v0.4.4
⌃ [f517fe37] Polyester v0.7.5
  [1d0040c9] PolyesterWeave v0.2.1
⌅ [f27b6e38] Polynomials v3.2.13
  [85a6dd25] PositiveFactorizations v0.2.4
  [d236fae5] PreallocationTools v0.4.12
⌃ [aea7be01] PrecompileTools v1.1.2
⌃ [21216c6a] Preferences v1.4.0
  [08abe8d2] PrettyTables v2.2.7
  [27ebfcd6] Primes v0.5.4
  [132c30aa] ProfileSVG v0.2.1
⌃ [92933f4c] ProgressMeter v1.7.2
⌃ [1fd47b50] QuadGK v2.8.2
⌃ [a08977f5] QuadratureRules v0.1.4
⌃ [c4ea9172] QuasiArrays v0.11.0
  [74087812] Random123 v1.6.1
  [fb686558] RandomExtensions v0.4.3
  [e6cf234a] RandomNumbers v1.5.3
  [c84ed2f1] Ratios v0.4.5
  [3cdcf5f2] RecipesBase v1.3.4
  [01d81517] RecipesPipeline v0.6.12
⌃ [731186ca] RecursiveArrayTools v2.38.7
⌃ [f2c3362d] RecursiveFactorization v0.2.18
  [189a3867] Reexport v1.2.2
  [05181044] RelocatableFolders v1.0.0
  [ae029012] Requires v1.3.0
  [ae5879a3] ResettableStacks v1.1.1
  [79098fc4] Rmath v0.7.1
  [47965b36] RootedTrees v2.19.2
⌅ [fb486d5c] RungeKutta v0.4.6
⌃ [7e49a35a] RuntimeGeneratedFunctions v0.5.11
  [fdea26ae] SIMD v3.4.5
  [94e857df] SIMDTypes v0.1.0
  [476501e8] SLEEFPirates v0.6.39
⌅ [0bca4576] SciMLBase v1.94.0
  [31c91b34] SciMLBenchmarks v0.1.3
⌃ [e9a6253c] SciMLNLSolve v0.1.8
  [c0aeaf25] SciMLOperators v0.3.6
  [6c6a2e73] Scratch v1.2.0
  [efcf1570] Setfield v1.1.1
  [992d4aef] Showoff v1.0.3
  [777ac1f9] SimpleBufferStream v1.1.0
⌃ [727e6d20] SimpleNonlinearSolve v0.1.19
⌅ [36b790f5] SimpleSolvers v0.2.4
  [699a6c99] SimpleTraits v0.9.4
  [ce78b400] SimpleUnPack v1.1.0
  [66db9d55] SnoopPrecompile v1.0.3
  [b85f4697] SoftGlobalScope v1.1.0
  [a2af1166] SortingAlgorithms v1.1.1
⌃ [47a9eef4] SparseDiffTools v2.4.1
  [e56a9233] Sparspak v0.3.9
⌃ [276daf66] SpecialFunctions v2.3.0
  [aedffcd0] Static v0.8.8
⌃ [0d7ed370] StaticArrayInterface v1.4.0
⌃ [90137ffa] StaticArrays v1.6.2
  [1e83bf80] StaticArraysCore v1.4.2
⌃ [82ae8749] StatsAPI v1.6.0
⌃ [2913bbd2] StatsBase v0.34.0
  [4c63d2b9] StatsFuns v1.3.0
  [7792a7ef] StrideArraysCore v0.4.17
  [69024149] StringEncodings v0.3.7
⌃ [892a3eda] StringManipulation v0.3.0
⌃ [c3572dad] Sundials v4.19.3
  [2efcf032] SymbolicIndexingInterface v0.2.2
⌃ [d1185830] SymbolicUtils v1.2.0
⌃ [0c5d862f] Symbolics v5.5.1
  [3783bdb8] TableTraits v1.0.1
⌃ [bd369af6] Tables v1.10.1
  [62fd8b95] TensorCore v0.1.1
  [8ea1fca8] TermInterface v0.3.3
  [8290d209] ThreadingUtilities v0.5.2
  [a759f4b9] TimerOutputs v0.5.23
  [c751599d] ToeplitzMatrices v0.8.2
  [0796e94c] Tokenize v0.5.25
  [3bb67fe8] TranscodingStreams v0.9.13
  [a2a6695c] TreeViews v0.3.0
  [d5829a12] TriangularSolve v0.1.19
  [410a4b4d] Tricks v0.1.7
  [781d530d] TruncatedStacktraces v1.4.0
⌃ [5c2747f8] URIs v1.4.2
  [3a884ed6] UnPack v1.0.2
  [1cfade01] UnicodeFun v0.4.1
⌃ [1986cc42] Unitful v1.16.0
  [45397f5d] UnitfulLatexify v1.6.3
  [a7c27f48] Unityper v0.1.5
  [41fe7b60] Unzip v0.2.0
  [3d5dd08c] VectorizationBase v0.21.64
  [81def892] VersionParsing v1.3.0
  [19fa3120] VertexSafeGraphs v0.2.0
  [44d3d7a6] Weave v0.10.12
  [efce3f68] WoodburyMatrices v0.5.5
  [ddb6d928] YAML v0.4.9
  [c2297ded] ZMQ v1.2.2
  [700de1a5] ZygoteRules v0.2.3
  [6e34b625] Bzip2_jll v1.0.8+0
  [83423d85] Cairo_jll v1.16.1+1
  [47200ebd] DecFP_jll v2.0.3+1
  [2e619515] Expat_jll v2.5.0+0
⌃ [b22a6f82] FFMPEG_jll v4.4.2+2
  [f5851436] FFTW_jll v3.3.10+0
  [34b6f7d7] FastTransforms_jll v0.6.2+0
  [a3f928ae] Fontconfig_jll v2.13.93+0
  [d7e528f0] FreeType2_jll v2.13.1+0
  [559328eb] FriBidi_jll v1.0.10+0
  [0656b61e] GLFW_jll v3.3.8+0
⌅ [d2c73de3] GR_jll v0.72.9+0
  [78b55507] Gettext_jll v0.21.0+0
  [f8c6e375] Git_jll v2.36.1+2
⌃ [7746bdde] Glib_jll v2.74.0+2
  [3b182d85] Graphite2_jll v1.3.14+0
⌃ [0234f1f7] HDF5_jll v1.12.2+2
  [2e76f6c2] HarfBuzz_jll v2.8.1+1
⌃ [1d5cc7b8] IntelOpenMP_jll v2023.1.0+0
  [aacddb02] JpegTurbo_jll v2.1.91+0
  [c1c5ebd0] LAME_jll v3.100.1+0
  [88015f11] LERC_jll v3.0.0+1
  [1d63c593] LLVMOpenMP_jll v15.0.4+0
  [aae0fff6] LSODA_jll v0.1.2+0
  [dd4b983a] LZO_jll v2.10.1+0
⌅ [e9f186c6] Libffi_jll v3.2.2+1
  [d4300ac3] Libgcrypt_jll v1.8.7+0
  [7e76a0d4] Libglvnd_jll v1.6.0+0
  [7add5ba3] Libgpg_error_jll v1.42.0+0
⌃ [94ce4f54] Libiconv_jll v1.16.1+2
  [4b2f31a3] Libmount_jll v2.35.0+0
  [89763e89] Libtiff_jll v4.5.1+1
  [38a345b3] Libuuid_jll v2.36.0+0
⌃ [856f044c] MKL_jll v2023.1.0+0
  [c771fb93] ODEInterface_jll v0.0.1+0
  [e7412a2a] Ogg_jll v1.3.5+1
⌅ [458c3c95] OpenSSL_jll v1.1.21+0
  [efe28fd5] OpenSpecFun_jll v0.5.5+0
  [91d4177d] Opus_jll v1.3.2+0
  [30392449] Pixman_jll v0.42.2+0
⌅ [c0090381] Qt6Base_jll v6.4.2+3
  [f50d1b31] Rmath_jll v0.4.0+0
⌅ [fb77eaff] Sundials_jll v5.2.1+0
⌃ [a2964d1f] Wayland_jll v1.21.0+0
  [2381bf8a] Wayland_protocols_jll v1.25.0+0
⌃ [02c8fc9c] XML2_jll v2.10.3+0
  [aed1982a] XSLT_jll v1.1.34+0
⌃ [ffd25f8a] XZ_jll v5.4.3+1
  [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.4+4
  [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.10+4
  [14d82f49] Xorg_libpthread_stubs_jll v0.1.1+0
  [c7cfdc94] Xorg_libxcb_jll v1.15.0+0
  [cc61e674] Xorg_libxkbfile_jll v1.1.2+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.4+0
  [3161d3a3] Zstd_jll v1.5.5+0
⌅ [214eeab7] fzf_jll v0.29.0+0
  [a4ae2306] libaom_jll v3.4.0+0
  [0ac62f75] libass_jll v0.15.1+0
  [f638f0a6] libfdk_aac_jll v2.0.2+0
  [b53b4c65] libpng_jll v1.6.38+0
  [a9144af2] libsodium_jll v1.0.20+0
  [f27f6e37] libvorbis_jll v1.3.7+1
  [1270edf5] x264_jll v2021.5.5+0
  [dfaa095f] x265_jll v3.5.0+0
⌃ [d8fb68d0] xkbcommon_jll v1.4.1+0
  [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.3
  [76f85450] LibGit2
  [8f399da3] Libdl
  [37e2e46d] LinearAlgebra
  [56ddb016] Logging
  [d6f4376e] Markdown
  [a63ad114] Mmap
  [ca575930] NetworkOptions v1.2.0
  [44cfe95a] Pkg v1.9.0
  [de0858da] Printf
  [9abbd945] Profile
  [3fa0cd96] REPL
  [9a3f8284] Random
  [ea8e919c] SHA v0.7.0
  [9e88b42a] Serialization
  [1a1011a3] SharedArrays
  [6462fe0b] Sockets
  [2f01184e] SparseArrays
  [10745b16] Statistics v1.9.0
  [4607b0f0] SuiteSparse
  [fa267f1f] TOML v1.0.3
  [a4e569a6] Tar v1.10.0
  [8dfed614] Test
  [cf7118a7] UUIDs
  [4ec0a83e] Unicode
  [e66e0078] CompilerSupportLibraries_jll v1.0.2+0
  [781609d7] GMP_jll v6.2.1+2
  [deac9b47] LibCURL_jll v7.84.0+0
  [29816b5a] LibSSH2_jll v1.10.2+0
  [3a97d323] MPFR_jll v4.1.1+4
  [c8ffd9c3] MbedTLS_jll v2.28.2+0
  [14a3606d] MozillaCACerts_jll v2022.10.11
  [4536629a] OpenBLAS_jll v0.3.21+4
  [05823500] OpenLibm_jll v0.8.1+0
  [efcefdf7] PCRE2_jll v10.42.0+0
  [bea87d4a] SuiteSparse_jll v5.10.1+6
  [83775a58] Zlib_jll v1.2.13+0
  [8e850b90] libblastrampoline_jll v5.8.0+0
  [8e850ede] nghttp2_jll v1.48.0+0
  [3f19e933] p7zip_jll v17.4.0+0
Info Packages marked with ⌃ and ⌅ have new versions available, 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.