diff --git a/docs/src/basics/InputOutput.md b/docs/src/basics/InputOutput.md index b1eb2905df..4568e10b0e 100644 --- a/docs/src/basics/InputOutput.md +++ b/docs/src/basics/InputOutput.md @@ -87,6 +87,77 @@ See [Symbolic Metadata](@ref symbolic_metadata). Metadata specified when creatin See [Linearization](@ref linearization). +## Real-time Input Handling During Simulation + +ModelingToolkit supports setting input values during simulation for variables marked with the `[input=true]` metadata. This is useful for real-time simulations, hardware-in-the-loop testing, interactive simulations, or any scenario where input values need to be determined during integration rather than specified beforehand. + +To use this functionality, variables must be marked as inputs using the `[input=true]` metadata and specified in the `inputs` keyword argument of `@mtkcompile`. + +There are two approaches to handling inputs during simulation: + +### Determinate Form: Using `Input` Objects + +When all input values are known beforehand, you can use the [`Input`](@ref) type to specify input values at specific time points. The solver will automatically apply these values using discrete callbacks. + +```@example inputs +using ModelingToolkit +using ModelingToolkit: t_nounits as t, D_nounits as D +using OrdinaryDiffEq +using Plots + +# Define system with an input variable +@variables x(t) [input=true] +@variables y(t) = 0 + +eqs = [D(y) ~ x] + +# Compile with inputs specified +@mtkcompile sys=System(eqs, t, [x, y], []) inputs=[x] + +prob = ODEProblem(sys, [], (0, 4)) + +# Create an Input object with predetermined values +input = Input(sys.x, [1, 2, 3, 4], [0, 1, 2, 3]) + +# Solve with the input - solver handles callbacks automatically +sol = solve(prob, [input], Tsit5()) + +plot(sol; idxs = [x, y]) +``` + +Multiple `Input` objects can be passed in a vector to handle multiple input variables simultaneously. + +### Indeterminate Form: Manual Input Setting with `set_input!` + +When input values need to be computed on-the-fly or depend on external data sources, you can manually set inputs during integration using [`set_input!`](@ref). This approach requires explicit control of the integration loop. + +```@example inputs +# Initialize the integrator +integrator = init(prob, Tsit5()) + +# Manually set inputs and step through time +set_input!(integrator, sys.x, 1.0) +step!(integrator, 1.0, true) + +set_input!(integrator, sys.x, 2.0) +step!(integrator, 1.0, true) + +set_input!(integrator, sys.x, 3.0) +step!(integrator, 1.0, true) + +set_input!(integrator, sys.x, 4.0) +step!(integrator, 1.0, true) + +# IMPORTANT: Must call finalize! to save all input callbacks +finalize!(integrator) + +plot(sol; idxs = [x, y]) +``` + +!!! warning "Always call `finalize!`" + + When using `set_input!`, you must call [`finalize!`](@ref) after integration is complete. This ensures that all discrete callbacks associated with input variables are properly saved in the solution. Without this call, input values may not be correctly recorded when querying the solution. + ## Docstrings ```@index @@ -96,4 +167,7 @@ Pages = ["InputOutput.md"] ```@docs; canonical=false ModelingToolkit.generate_control_function ModelingToolkit.build_explicit_observed_function +ModelingToolkit.Input +ModelingToolkit.set_input! +ModelingToolkit.finalize! ``` diff --git a/src/ModelingToolkit.jl b/src/ModelingToolkit.jl index 6c79eb4fb1..d39d657390 100644 --- a/src/ModelingToolkit.jl +++ b/src/ModelingToolkit.jl @@ -165,6 +165,7 @@ include("constants.jl") include("utils.jl") +export set_input!, finalize!, Input include("systems/index_cache.jl") include("systems/parameter_buffer.jl") include("systems/abstractsystem.jl") @@ -175,6 +176,7 @@ include("systems/state_machines.jl") include("systems/analysis_points.jl") include("systems/imperative_affect.jl") include("systems/callbacks.jl") +include("systems/inputs.jl") include("systems/system.jl") include("systems/codegen_utils.jl") include("problems/docs.jl") diff --git a/src/systems/abstractsystem.jl b/src/systems/abstractsystem.jl index 51a495e4ff..1c00079a9f 100644 --- a/src/systems/abstractsystem.jl +++ b/src/systems/abstractsystem.jl @@ -789,7 +789,8 @@ const SYS_PROPS = [:eqs :index_cache :isscheduled :costs - :consolidate] + :consolidate + :input_functions] for prop in SYS_PROPS fname_get = Symbol(:get_, prop) diff --git a/src/systems/inputs.jl b/src/systems/inputs.jl new file mode 100644 index 0000000000..413ab5dbf7 --- /dev/null +++ b/src/systems/inputs.jl @@ -0,0 +1,196 @@ +""" + Input(var, data::Vector{<:Real}, time::Vector{<:Real}) + +Create an `Input` object that specifies predetermined input values for a variable at specific time points. + +# Arguments +- `var`: The symbolic variable (marked with `[input=true]` metadata) to be used as an input. +- `data`: A vector of real values that the input variable should take at the corresponding time points. +- `time`: A vector of time points at which the input values should be applied. Must be the same length as `data`. + +# Description +The `Input` struct is used with the extended `solve` method to provide time-varying inputs to a system +during simulation. When passed to `solve(prob, [input1, input2, ...], alg)`, the solver will automatically +set the input variable to the specified values at the specified times using discrete callbacks. + +This provides a "determinate form" of input handling where all input values are known a priori, +as opposed to setting inputs manually during integration with [`set_input!`](@ref). + +See also [`set_input!`](@ref), [`finalize!`](@ref) +""" +struct Input + var::Num + data::SVector + time::SVector +end + +function Input(var, data::Vector{<:Real}, time::Vector{<:Real}) + n = length(data) + return Input(var, SVector{n}(data), SVector{n}(time)) +end + +struct InputFunctions{S, O} + events::Tuple{SymbolicDiscreteCallback} + vars::Tuple{SymbolicUtils.BasicSymbolic{Real}} + setters::Tuple{SymbolicIndexingInterface.ParameterHookWrapper{S, O}} +end + +function InputFunctions(events::Vector, vars::Vector, setters::Vector) + InputFunctions(Tuple(events), Tuple(vars), Tuple(setters)) +end + +""" + set_input!(integrator, var, value::Real) + +Set the value of an input variable during integration. + +# Arguments +- `integrator`: An ODE integrator object (from `init(prob, alg)` or available in callbacks). +- `var`: The symbolic input variable to set (must be marked with `[input=true]` metadata and included in the `inputs` keyword of `@mtkcompile`). +- `value`: The new real-valued input to assign to the variable. +- `input_funs` (optional): The `InputFunctions` object associated with the system. If not provided, it will be retrieved from `integrator.f.sys`. + +# Description +This function allows you to manually set input values during integration, providing an "indeterminate form" +of input handling where inputs can be computed on-the-fly. This is useful when input values depend on +runtime conditions, external data sources, or interactive user input. + +After setting input values with `set_input!`, you must call [`finalize!`](@ref) at the end of integration +to ensure all discrete callbacks are properly saved. + +# Example +```julia +@variables x(t) [input=true] +@variables y(t) = 0 + +eqs = [D(y) ~ x] +@mtkcompile sys = System(eqs, t, [x, y], []) inputs=[x] + +prob = ODEProblem(sys, [], (0, 4)) +integrator = init(prob, Tsit5()) + +# Set input and step forward +set_input!(integrator, sys.x, 1.0) +step!(integrator, 1.0, true) + +set_input!(integrator, sys.x, 2.0) +step!(integrator, 1.0, true) + +# Must call finalize! at the end +finalize!(integrator) +``` + +See also [`finalize!`](@ref), [`Input`](@ref) +""" +function set_input!(input_funs::InputFunctions, integrator::OrdinaryDiffEqCore.ODEIntegrator, var, value::Real) + i = findfirst(isequal(var), input_funs.vars) + setter = input_funs.setters[i] + event = input_funs.events[i] + + setter(integrator, value) + save_callback_discretes!(integrator, event) + u_modified!(integrator, true) + return nothing +end +function set_input!(integrator, var, value::Real) + set_input!(get_input_functions(integrator.f.sys), integrator, var, value) +end + +""" + finalize!(integrator) + +Finalize all input callbacks at the end of integration. + +# Arguments +- `integrator`: An ODE integrator object (from `init(prob, alg)` or available in callbacks). +- `input_funs` (optional): The `InputFunctions` object associated with the system. If not provided, it will be retrieved from `integrator.f.sys`. + +# Description +This function must be called after using [`set_input!`](@ref) to manually set input values during integration. +It ensures that all discrete callbacks associated with input variables are properly saved in the solution, +making the input values accessible when querying the solution at specific time points. + +Without calling `finalize!`, input values set with `set_input!` may not be correctly recorded in the +final solution object, leading to incorrect results when indexing the solution. + +See also [`set_input!`](@ref), [`Input`](@ref) +""" +function finalize!(input_funs::InputFunctions, integrator) + for i in eachindex(input_funs.vars) + save_callback_discretes!(integrator, input_funs.events[i]) + end + + return nothing +end +finalize!(integrator) = finalize!(get_input_functions(integrator.f.sys), integrator) + +function (input_funs::InputFunctions)(integrator, var, value::Real) + set_input!(input_funs, integrator, var, value) +end +(input_funs::InputFunctions)(integrator) = finalize!(input_funs, integrator) + +function build_input_functions(sys, inputs) + + # Here we ensure the inputs have metadata marking the discrete variables as parameters. In some + # cases the inputs can be fed to this function before they are converted to parameters by mtkcompile. + vars = SymbolicUtils.BasicSymbolic[isparameter(x) ? x : toparam(x) + for x in unwrap.(inputs)] + setters = [] + events = SymbolicDiscreteCallback[] + defaults = get_defaults(sys) + if !isempty(vars) + for x in vars + affect = ImperativeAffect((m, o, c, i)->m, modified = (; x)) + sdc = SymbolicDiscreteCallback(Inf, affect) + + push!(events, sdc) + + # ensure that the ODEProblem does not complain about missing parameter map + if !haskey(defaults, x) + push!(defaults, x => 0.0) + end + end + + @set! sys.discrete_events = events + @set! sys.index_cache = ModelingToolkit.IndexCache(sys) + @set! sys.defaults = defaults + + setters = [SymbolicIndexingInterface.setsym(sys, x) for x in vars] + + @set! sys.input_functions = InputFunctions(events, vars, setters) + end + + return sys +end + +function CommonSolve.solve(prob::SciMLBase.AbstractDEProblem, inputs::Vector{Input}, args...; kwargs...) + tstops = Float64[] + callbacks = DiscreteCallback[] + + # set_input! + for input::Input in inputs + tstops = union(tstops, input.time) + condition = (u, t, integrator) -> any(t .== input.time) + affect! = function (integrator) + @inbounds begin + i = findfirst(integrator.t .== input.time) + set_input!(integrator, input.var, input.data[i]) + end + end + push!(callbacks, DiscreteCallback(condition, affect!)) + + # DiscreteCallback doesn't hit on t==0, workaround... + if input.time[1] == 0 + prob.ps[input.var] = input.data[1] + end + end + + # finalize! + t_end = prob.tspan[2] + condition = (u, t, integrator) -> (t == t_end) + affect! = (integrator) -> finalize!(integrator) + push!(callbacks, DiscreteCallback(condition, affect!)) + push!(tstops, t_end) + + return solve(prob, args...; tstops, callback = CallbackSet(callbacks...), kwargs...) +end diff --git a/src/systems/system.jl b/src/systems/system.jl index 12f5d1d250..87f81403cd 100644 --- a/src/systems/system.jl +++ b/src/systems/system.jl @@ -259,6 +259,11 @@ struct System <: IntermediateDeprecationSystem The `Schedule` containing additional information about the simplified system. """ schedule::Union{Schedule, Nothing} + """ + $INTERNAL_FIELD_WARNING + Functions used to set input variables with `set_input!` and `finalize!` functions + """ + input_functions::Union{InputFunctions, Nothing} function System( tag, eqs, noise_eqs, jumps, constraints, costs, consolidate, unknowns, ps, @@ -271,7 +276,7 @@ struct System <: IntermediateDeprecationSystem complete = false, index_cache = nothing, ignored_connections = nothing, preface = nothing, parent = nothing, initializesystem = nothing, is_initializesystem = false, is_discrete = false, isscheduled = false, - schedule = nothing; checks::Union{Bool, Int} = true) + schedule = nothing, input_functions = nothing; checks::Union{Bool, Int} = true) if is_initializesystem && iv !== nothing throw(ArgumentError(""" Expected initialization system to be time-independent. Found independent @@ -310,7 +315,7 @@ struct System <: IntermediateDeprecationSystem tstops, inputs, outputs, tearing_state, namespacing, complete, index_cache, ignored_connections, preface, parent, initializesystem, is_initializesystem, is_discrete, - isscheduled, schedule) + isscheduled, schedule, input_functions) end end diff --git a/src/systems/systems.jl b/src/systems/systems.jl index 49efb3e1d4..a5b0325129 100644 --- a/src/systems/systems.jl +++ b/src/systems/systems.jl @@ -51,6 +51,11 @@ function mtkcompile( @set! newsys.parent = complete(sys; split = false, flatten = false) end newsys = complete(newsys; split) + + if !isempty(inputs) + newsys = build_input_functions(newsys, inputs) + end + if newsys′ isa Tuple idxs = [parameter_index(newsys, i) for i in io[1]] return newsys, idxs diff --git a/test/inputs.jl b/test/inputs.jl new file mode 100644 index 0000000000..eb9314edd2 --- /dev/null +++ b/test/inputs.jl @@ -0,0 +1,59 @@ +using ModelingToolkit +using ModelingToolkit: t_nounits as t, D_nounits as D +using OrdinaryDiffEq +using Test + +# ----------------------------------------- +# ----- example --------------------------- +# ----------------------------------------- + +vars = @variables begin + x(t), [input=true] + + # states + y(t) = 0 +end + +eqs = [ + D(y) ~ x +] + +@mtkcompile sys=System(eqs, t, vars, []) inputs=[x] + +@test !isnothing(ModelingToolkit.get_input_functions(sys)) + +prob = ODEProblem(sys, [], (0, 4)) + +# indeterminate form ----------------------- + +integrator = init(prob, Tsit5()) + +set_input!(integrator, sys.x, 1.0) +step!(integrator, 1.0, true) + +set_input!(integrator, sys.x, 2.0) +step!(integrator, 1.0, true) + +set_input!(integrator, sys.x, 3.0) +step!(integrator, 1.0, true) + +set_input!(integrator, sys.x, 4.0) +step!(integrator, 1.0, true) + +finalize!(integrator) + +@test integrator.sol(0.0; idxs = sys.x) == 1.0 +@test integrator.sol(1.0; idxs = sys.x) == 2.0 +@test integrator.sol(2.0; idxs = sys.x) == 3.0 +@test integrator.sol(3.0; idxs = sys.x) == 4.0 +@test integrator.sol(4.0; idxs = sys.y) ≈ 10.0 + +# determinate form ----------------------- +input = Input(sys.x, [1, 2, 3, 4], [0, 1, 2, 3]) +sol = solve(prob, [input], Tsit5()); + +@test sol(0.0; idxs = sys.x) == 1.0 +@test sol(1.0; idxs = sys.x) == 2.0 +@test sol(2.0; idxs = sys.x) == 3.0 +@test sol(3.0; idxs = sys.x) == 4.0 +@test sol(4.0; idxs = sys.y) ≈ 10.0 diff --git a/test/runtests.jl b/test/runtests.jl index 522470b896..344ee2f595 100644 --- a/test/runtests.jl +++ b/test/runtests.jl @@ -34,6 +34,7 @@ end @safetestset "Direct Usage Test" include("direct.jl") @safetestset "System Linearity Test" include("linearity.jl") @safetestset "Input Output Test" include("input_output_handling.jl") + @safetestset "Inputs" include("inputs.jl") @safetestset "Clock Test" include("clock.jl") @safetestset "ODESystem Test" include("odesystem.jl") @safetestset "Dynamic Quantities Test" include("dq_units.jl") @@ -103,7 +104,7 @@ end @safetestset "Fractional Differential Equations Tests" include("fractional_to_ordinary.jl") end end - + if GROUP == "All" || GROUP == "SymbolicIndexingInterface" @safetestset "SymbolicIndexingInterface test" include("symbolic_indexing_interface.jl") @safetestset "SciML Problem Input Test" include("sciml_problem_inputs.jl")