Skip to content

Large allocations with MTK event #4015

@klinders

Description

@klinders

Describe the bug 🐞
Adding an event to the MTK model will slow down the solve significantly due to excessive allocations.

Expected behavior

Minimal performance impact when using events.

Minimal Reproducible Example 👇

using ModelingToolkit
using ModelingToolkit: t_nounits as t, D_nounits as D

function affect!(mod,obs,ctx,int)
    ModelingToolkit.terminate!(int)

    return (;)
end

function FOL_event(;name)
    @parameters begin
        τ = 3.0 # parameters
    end
    @variables begin
        x(t) = 0.0 # dependent variables
    end
    eq=[
        D(x) ~ (1 - x) / τ
    ]

    events = [x ~ 0.6]=>(affect!,(;))

    return System(eq, t; name=name, continuous_events=events)
end

function FOL(;name)
    @parameters begin
        τ = 3.0 # parameters
    end
    @variables begin
        x(t) = 0.0 # dependent variables
    end
    eq=[
        D(x) ~ (1 - x) / τ
    ]

    return System(eq, t; name=name)
end

using OrdinaryDiffEq
@mtkbuild fol = FOL()
@mtkbuild fol_event = FOL_event()

print("Starting problems\n")
print("No event:")
@time prob = ODEProblem(fol, [], (0.0, 10.0), dense=true,save_everystep=false)
print("With event:")
@time prob_event = ODEProblem(fol_event, [], (0.0, 10.0), dense=true, save_everystep=false)

print("Solving problems\n")
print("No event:")
@time sol = solve(prob,Tsit5())
print("With event:")
@time sol_event = solve(prob_event, Tsit5())

Simulation output ⚠️
I ran the code above in a REPL twice to get rid of the compilation time.

Starting problems
No event:  0.009388 seconds (39.02 k allocations: 2.201 MiB)
With event:  0.007077 seconds (36.27 k allocations: 1.683 MiB)
Solving problems
No event:  0.000044 seconds (88 allocations: 6.578 KiB)
With event:  0.119570 seconds (144.21 k allocations: 8.894 MiB, 99.95% compilation time: 100% of which was recompilation)

Environment (please complete the following information):

  • Output of using Pkg; Pkg.status()
Status `~/Documents/Simulations/Julia/BPackSimulator/Project.toml`
  [e9eeebe6] BatteryPeck v0.1.0 `src/BatteryPeck`
  [634d3b9d] DrWatson v2.19.1
  [961ee093] ModelingToolkit v10.26.1
  [16a59e39] ModelingToolkitStandardLibrary v2.25.0
  [1dea7af3] OrdinaryDiffEq v6.103.0
  [91a5bcdd] Plots v1.41.1
  [295af30f] Revise v3.12.0
  • Output of versioninfo()
Julia Version 1.10.10
Commit 95f30e51f41 (2025-06-27 09:51 UTC)
Build Info:
  Official https://julialang.org/ release
Platform Info:
  OS: Linux (x86_64-linux-gnu)
  CPU: 12 × 13th Gen Intel(R) Core(TM) i7-1365U
  WORD_SIZE: 64
  LIBM: libopenlibm
  LLVM: libLLVM-15.0.7 (ORCJIT, goldmont)
Threads: 1 default, 0 interactive, 1 GC (on 12 virtual cores)

Additional context
The topic seems similar to Significant allocations with Callbacks (Tsit5) , but since those fixes were merged a long time ago, I decided to open a new issue.

Metadata

Metadata

Assignees

No one assigned

    Labels

    bugSomething isn't working

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions