Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Errors in affect codegen #2994

Open
baggepinnen opened this issue Aug 28, 2024 · 14 comments
Open

Errors in affect codegen #2994

baggepinnen opened this issue Aug 28, 2024 · 14 comments
Labels
bug Something isn't working events

Comments

@baggepinnen
Copy link
Contributor

The following system simulates correctly without the event. The event affect function of the event is empty, so I would have expected it to have no effect.

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

@component function ContactForce2(; name, surface=nothing)
    vars = @variables begin
        q(t) = 1
        v(t) = 0
        f(t)
        h(t)
        hd(t)
        hdd(t)
        (λ(t)=0), [irreducible=true]
    end
    pars = @parameters begin
        contact::Int = 0    # discrete.time state variable
        a0 = 100
        a1 = 20
        a2 = 1e6
    end

    equations = [
        0 ~ ifelse((contact == 1), hdd + a1*hd + a0*h + a2*h^3, λ)
        f ~ contact*λ
        D(q) ~ v

        1 * D(v) ~ -1 * 9.81 + f

        h ~ q
        hd ~ D(h)
        hdd ~ D(hd)
    ]

    function affect!(integ, u, p, _)
    end
    continuous_events = [h ~ 0] => (affect!, [h], [], [], nothing)

    ODESystem(equations, t, vars, pars; name, continuous_events)
end
@named model = ContactForce2()
model = complete(model)
ssys = structural_simplify(model)
prob = ODEProblem(ssys, [], (0.0, 5.0))
sol = solve(prob, Rodas4())

using Plots
plot(sol, layout=5)

without event an object is falling freely due to gravity
image

with empty event
image

@baggepinnen baggepinnen added bug Something isn't working events labels Aug 28, 2024
@baggepinnen
Copy link
Contributor Author

The problem appears related to the equation

0 ~ ifelse((contact == 1), hdd + a1*hd + a0*h + a2*h^3, λ)

Since contact === 0 in this simulation, this equation always reads 0 ~ λ. If I use 0 ~ λ as equation instead, it works as expected:

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

@component function ContactForce2(; name, surface=nothing)
    vars = @variables begin
        q(t) = 1
        v(t) = 0
        f(t)
        (h(t) = 1), [irreducible=true]
        hd(t)
        hdd(t)
        (λ(t)=0), [irreducible=true]
    end
    pars = @parameters begin
        contact::Int = 0    # discrete.time state variable
        a0 = 100
        a1 = 20
        a2 = 1e6
    end

    equations = [
        0 ~ λ# ifelse((contact == 1), hdd + a1*hd + a0*h + a2*h^3, λ)
        f ~ contact*λ
        D(q) ~ v

        1 * D(v) ~ -1 * 9.81 + f

        h ~ q
        hd ~ D(h)
        hdd ~ D(hd)
    ]

    function affect!(integ, u, p, _)
    end
    continuous_events = [h ~ 0] => (affect!, [h], [], [], nothing)

    ODESystem(equations, t, vars, pars; name, continuous_events)
end
@named model = ContactForce2()
model = complete(model)
ssys = structural_simplify(model)
prob = ODEProblem(ssys, [], (0.0, 5.0))
sol = solve(prob, Rodas4())

using Plots
plot(sol, layout=5)

image

@baggepinnen
Copy link
Contributor Author

baggepinnen commented Aug 28, 2024

Even more weird stuff going on. If I set dtmax to something small, the h coordinate is reset to it's initial value each time the empty event triggers

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

@component function ContactForce2(; name, surface=nothing)
    vars = @variables begin
        (q(t) = 1), [irreducible=true]
        v(t) = 0
        (f(t)=0), [irreducible=true]
        (h(t) = 1), [irreducible=true]
        (hd(t) = 0)#, [irreducible=true]
        hdd(t)#, [irreducible=true]
        (λ(t)=0), [irreducible=true]
    end
    pars = @parameters begin
        contact::Int = 0    # discrete.time state variable
        a0 = 100
        a1 = 20
        a2 = 0*1e6
    end

    equations = [
        0 ~ ifelse((contact == 1), hdd + a1*hd + a0*h + a2*h^3, λ)
        f ~ contact*λ
        D(q) ~ v

        1 * D(v) ~ -1 * 9.81 + f

        h ~ q
        hd ~ D(h)
        hdd ~ D(hd)
    ]

    function affect!(integ, u, p, _)
    end
    continuous_events = [h ~ 0] => (affect!, [h], [], [], nothing)

    ODESystem(equations, t, vars, pars; name, continuous_events)
end


@named model = ContactForce2()
model = complete(model)
ssys = structural_simplify(model)
prob = ODEProblem(ssys, [], (0.0, 2.0))
sol = solve(prob, Rodas4(), dtmax=0.001)

using Plots
plot(sol, layout=6, size=(1000, 1000))

image

@baggepinnen baggepinnen changed the title Presence of empty event screws up result Errors in affect codegen Aug 28, 2024
@ChrisRackauckas
Copy link
Member

@BenChung can you take this one?

@BenChung
Copy link
Contributor

This is the bug we talked about on Slack where we reinitialize the system after the affect fires using the initialization system for the initial condition, rather than the "running" condition. I probably can take it - eventually - but I don't know enough about the initialization system to really be able to fix it properly.

@ChrisRackauckas ChrisRackauckas assigned BenChung and unassigned BenChung Aug 28, 2024
@ChrisRackauckas
Copy link
Member

Oh I see. This needs to be handled through what I had mentioned with the ImplicitDiscreteSystem thing if we're to make it robust.

@BenChung
Copy link
Contributor

Yeah. My feeling - for now - is that it's better to not run initialization (or error if there is an initialization) after an affect so that the user doesn't get surprised by it and takes on the responsibility for ensuring that the requirements of a DAE are satisfied, and then we need to revisit this problem once we have better handling for system "live" initialization.

@ChrisRackauckas
Copy link
Member

That's simply impossible.

@baggepinnen
Copy link
Contributor Author

Why is running the initialization system required at all? If there are no algebraic equations, the state cannot be in feasible. If there are algebraic equations, the DAE solver will have to find a root to progress anyways and will solve the problem then?

The initialization system would be the incorrect system to solve anyways, parameters might have changed and equations for initial condition should no longer be included, it would be a completely different system.

@ChrisRackauckas
Copy link
Member

You always have to run initialization after any affect! because otherwise you cannot guarantee a consistent state, and if you don't have a consistent state then you cannot take a step of the integrator.

The initialization system would be the incorrect system to solve anyways, parameters might have changed and equations for initial condition should no longer be included, it would be a completely different system.

That's the point of the implicit discrete form.

@HKruenaegel
Copy link

I have the same probleme with a discrete Callback, states are set to initial value every time the callback is fired.

What can I do to get the same behaviour as with MTK v8.73?

@HKruenaegel
Copy link

I have a model which is developed with MTK v8.72.0. With v8.75.0 it becomes slow an interrupts without in error after 3 seconds. With v9.19.0 it becomes much slower and interrupts with error at the same time. With >v9.36 it has the Problem described above.

I will try to shrink down the system for a example that I can share. But it seems that it's something with the descrete callback. Are there similar known Issues yet?

@BenChung
Copy link
Contributor

@HKruenaegel We should be fixing the incorrect initialization system issue that causes the system to go to all-0s after an affect fires soon(^tm) when my development branch gets merged. I haven't had an opportunity to look into performance, though, so can't comment on that.

@HKruenaegel
Copy link

I tried again to run my model with MTKv9.30.0. If I give all states in u0 and the guesses are empty the model works, but as mentioned before it is 30 times slower than with MTKv8.73.2.

Is that a general observation with MTK9, @ChrisRackauckas?

@HKruenaegel
Copy link

Here is another example for the problem initially mentioned in this issue. This is a system which represents this circuit:
image

The resistor Rd is a ideal diode switch, therfore a Callback is defined, which sets the resistance to 1e4 if the voltage over Rd is <0 and to 1e-3 if it is >0.

The result of the simulation gives:
image

The current through the inductor is correct, also the voltage of the capacitor v3 during conducting state of the diode. But when the callback for switching off the diode is fired the capacitor voltage v3 is set to zero, which is not correct.

For this example MTK v.9.40.0 is used.

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

@parameters Rd Rsw C1 L1 V1 Rl
@variables I_L1(t) I_V1(t) v1(t) v2(t) [irreducible = true] v3(t) [irreducible = true]

eqs =         [I_L1 + I_V1~ 0
-I_L1 + v2/Rsw + v2/Rd - v3/Rd~ 0
C1*D(v3) + v3/Rl - v2/Rd + v3/Rd~ 0
v1~ 10*sin(2*pi*50*t)
-D(I_L1)*L1 + v1 - v2~ 0]

function D_on(i, u, p,ctx)
    i.ps[p.Rd]=1e-3
end
function D_off(i, u, p,ctx)
    i.ps[p.Rd]=1e3
end
c = ModelingToolkit.SymbolicContinuousCallback(
    [v2-v3~0], (D_on, [v2,v3], [Rd], [], nothing);affect_neg =(D_off, [v2,v3], [Rd], [], nothing),
    rootfind = SciMLBase.LeftRootFind)

@mtkbuild pend = ODESystem(eqs, t;continuous_events=c)

u0 = [I_L1=>0.0,
v3=>0.0,
]

g=[v1=>0.0,
v2=>0.0,
I_V1=>0.0]
p = [        Rd => 0.001,
Rsw=>1e6,
C1=>0.01,
L1=>0.001,
V1=>10.0,
Rl=>20.0]
prob = ODEProblem(pend, u0, (0.0, 1.5), p, guesses = g)

@time sol = solve(prob, Rodas5P(),dtmax=1e-4)

Sign up for free to join this conversation on GitHub. Already have an account? Sign in to comment
Labels
bug Something isn't working events
Projects
None yet
Development

No branches or pull requests

5 participants
@BenChung @ChrisRackauckas @baggepinnen @HKruenaegel and others