Skip to content

Conversation

@ChrisRackauckas
Copy link
Member

Fixes https://discourse.julialang.org/t/function-compiles-every-time/63273/5 by giving the arguments unique names. MWE:

using ModelingToolkit, OrdinaryDiffEq

@variables t x(t) RHS(t)  # independent and dependent variables
@parameters τ       # parameters
D = Differential(t) # define an operator for the differentiation w.r.t. time

# your first ODE, consisting of a single equation, indicated by ~
@named fol_separate = ODESystem([ RHS  ~ (1 - x)/τ,
                                  D(x) ~ RHS ])

sys = structural_simplify(fol_separate)

@time begin
    prob = ODEProblem(sys, [x => 0.0], (0.0,10.0), [τ => 3.0])
    sol = solve(prob, Tsit5())
end
  0.001047 seconds (2.20 k allocations: 126.812 KiB)
  0.001023 seconds (2.20 k allocations: 126.812 KiB)
  0.001021 seconds (2.20 k allocations: 125.672 KiB)
  0.001045 seconds (2.20 k allocations: 125.672 KiB)

while it used to be 99.9% compile time. Tested on <1 and rebased.

Fixes https://discourse.julialang.org/t/function-compiles-every-time/63273/5 by giving the arguments unique names. MWE:

```julia
using ModelingToolkit, OrdinaryDiffEq

@variables t x(t) RHS(t)  # independent and dependent variables
@parameters τ       # parameters
D = Differential(t) # define an operator for the differentiation w.r.t. time

# your first ODE, consisting of a single equation, indicated by ~
@nAmed fol_separate = ODESystem([ RHS  ~ (1 - x)/τ,
                                  D(x) ~ RHS ])

sys = structural_simplify(fol_separate)

@time begin
    prob = ODEProblem(sys, [x => 0.0], (0.0,10.0), [τ => 3.0])
    sol = solve(prob, Tsit5())
end
```

```
  0.001047 seconds (2.20 k allocations: 126.812 KiB)
  0.001023 seconds (2.20 k allocations: 126.812 KiB)
  0.001021 seconds (2.20 k allocations: 125.672 KiB)
  0.001045 seconds (2.20 k allocations: 125.672 KiB)
```

while it used to be 99.9% compile time.
linenumbers = true)

dargs = map(arg -> destructure_arg(arg, !checkbounds), [args...])
dargs = map((x) -> destructure_arg(x[2], !checkbounds, Symbol("ˍ₋arg$(x[1])")), enumerate([args...]))
Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

two unicode underscores are better than one?

@codecov-commenter
Copy link

codecov-commenter commented Jun 22, 2021

⚠️ Please install the 'codecov app svg image' to ensure uploads and comments are reliably processed by Codecov.

Codecov Report

All modified and coverable lines are covered by tests ✅

Project coverage is 74.36%. Comparing base (1d89d76) to head (e9810a4).
Report is 1880 commits behind head on master.

❗ Your organization needs to install the Codecov GitHub app to enable full functionality.

Additional details and impacted files
@@            Coverage Diff             @@
##           master     #272      +/-   ##
==========================================
- Coverage   74.44%   74.36%   -0.09%     
==========================================
  Files          19       19              
  Lines        1761     1759       -2     
==========================================
- Hits         1311     1308       -3     
- Misses        450      451       +1     

☔ View full report in Codecov by Sentry.
📢 Have feedback on the report? Share it here.


function numbered_expr(de::Equation,varnumbercache,args...;varordering = args[1],
lhsname=gensym("du"),rhsnames=[gensym("MTK") for i in 1:length(args)],offset=0)
lhsname=:du,rhsnames=[Symbol("MTK$i") for i in 1:length(args)],offset=0)
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

These are only used in C and MATLAB outputs iirc

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

Might as well make the codegen stable there too.

Copy link
Member Author

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

And in that case, we really don't need to gensym since it's not going to take globals from Julia 😅

function destructure_arg(arg::Union{AbstractArray, Tuple}, inbounds, name)
if !(arg isa Arr)
DestructuredArgs(map(value, arg), inbounds=inbounds)
DestructuredArgs(map(value, arg), name, inbounds=inbounds)
Copy link
Member

Choose a reason for hiding this comment

The reason will be displayed to describe this comment to others. Learn more.

nice!

@YingboMa
Copy link
Member

This might not solve the root issue. MTK is still regenerating the same function again.

@ChrisRackauckas
Copy link
Member Author

That function generation is rather cheap though. This is orders of magnitude better than before at least.

@YingboMa YingboMa merged commit 584721f into master Jun 22, 2021
@YingboMa YingboMa deleted the gensym branch June 22, 2021 18:21
Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment

Labels

None yet

Projects

None yet

Development

Successfully merging this pull request may close these issues.

5 participants