Skip to content

Discrete-time variables not found when forming InitializationProblem #2699

Open
@baggepinnen

Description

@baggepinnen

The following example fails with

julia> prob = ODEProblem(m, [m.delay.u(k-3)=>0, m.delay.u(k-2)=>0, m.delay.u(k-1)=>0], (0.0, 10.0))
ERROR: KeyError: key Shift(nothing, -3)(delay₊u(t)) not found
Stacktrace:
  [1] getindex(h::Dict{Any, Any}, key::SymbolicUtils.BasicSymbolic{Real})
    @ Base ./dict.jl:498
  [2] TearingState(sys::NonlinearSystem; quick_cancel::Bool, check::Bool)
    @ ModelingToolkit ~/.julia/dev/ModelingToolkit/src/systems/systemstructure.jl:393
  [3] TearingState(sys::NonlinearSystem)
    @ ModelingToolkit ~/.julia/dev/ModelingToolkit/src/systems/systemstructure.jl:250
  [4] __structural_simplify(sys::NonlinearSystem, io::Nothing; simplify::Bool, kwargs::@Kwargs{fully_determined::Bool})
    @ ModelingToolkit ~/.julia/dev/ModelingToolkit/src/systems/systems.jl:58
  [5] __structural_simplify
    @ ~/.julia/dev/ModelingToolkit/src/systems/systems.jl:55 [inlined]
  [6] structural_simplify(sys::NonlinearSystem, io::Nothing; simplify::Bool, split::Bool, kwargs::@Kwargs{})
    @ ModelingToolkit ~/.julia/dev/ModelingToolkit/src/systems/systems.jl:22
  [7] structural_simplify (repeats 2 times)
    @ ~/.julia/dev/ModelingToolkit/src/systems/systems.jl:19 [inlined]
  [8] ModelingToolkit.InitializationProblem{…}(sys::ODESystem, t::Float64, u0map::Vector{…}, parammap::SciMLBase.NullParameters; guesses::Dict{…}, check_length::Bool, warn_initialize_determined::Bool, kwargs::@Kwargs{})

despite m.delay.u(k-3)=>0 being defined.

using ModelingToolkit
using ModelingToolkitStandardLibrary.Blocks
k = ShiftIndex()

@mtkmodel Del begin
    @extend u, y = siso = SISO()
    @structural_parameters begin
        n = 1
    end
    @equations begin
        y ~ u(k-n)
    end
end
@mtkmodel DiscreteInt begin
    @extend u, y = siso = SISO()
    @variables begin
        x(t) = 0.0, [description = "State of Integrator"]
    end
    @equations begin
        x(k) ~ x(k-1) + u(k-1)
        y ~ x
    end
end
@mtkmodel DelayModel begin
    @components begin
        delay = Del(n = 3)
        plant = DiscreteInt(x=0)
        input = Constant(k = 1)
    end
    @equations begin
        connect(input.output, delay.input)
        connect(delay.output, plant.input)
    end
end

@mtkbuild m = DelayModel()
prob = ODEProblem(m, [m.delay.u(k-3)=>0, m.delay.u(k-2)=>0, m.delay.u(k-1)=>0], (0.0, 10.0))

If the example is simplified to

@mtkmodel Del begin
    @extend u, y = siso = SISO()
    @structural_parameters begin
        n = 1
    end
    @equations begin
        y ~ u(k-n)
    end
end

@mtkmodel DelayModel begin
    @components begin
        delay = Del(n = 3)
        input = Constant(k = 1)
    end
    @equations begin
        connect(input.output, delay.input)
    end
end

@mtkbuild m = DelayModel()
prob = ODEProblem(m, [m.delay.u(k-3)=>0, m.delay.u(k-2)=>0, m.delay.u(k-1)=>0], (0.0, 10.0))

The initialization does work correctly, even though the simplification removed the system plant = DiscreteInt() rather than the system delay

Metadata

Metadata

Assignees

No one assigned

    Labels

    Type

    No type

    Projects

    No projects

    Milestone

    No milestone

    Relationships

    None yet

    Development

    No branches or pull requests

    Issue actions