FSimBase.jl is the lightweight base library for numerical simulation supporting nested dynamical systems and macro-based data logger, compatible with DifferentialEquations.jl.
- FSimBase.jl works alone! For more functionality, see FlightSims.jl.
- In FSimBase.jl,
you must specify the differential equation problem and the differential equation solver.
Only
ODEProblem
andDiscreteProblem
are tested.
Main APIs are provided in src/APIs
.
AbstractEnv
: an abstract type for user-defined and predefined environments. In general, environments is a sub-type ofAbstractEnv
.struct LinearSystemEnv <: AbstractEnv A B end
State(env::AbstractEnv)
: return a function that produces structured states.function State(env::LinearSystemEnv) @unpack B = env n = size(B)[1] return function (x) @assert length(x) == n x end end
Dynamics!(env::AbstractEnv)
: return a function that maps in-place dynamics, compatible with DifferentialEquations.jl. User can extend these methods or simply define other methods.function Dynamics!(env::LinearSystemEnv) @unpack A, B = env @Loggable function dynamics!(dx, x, p, t; u) # data would not be saved without @Loggable. Follow this form! @log state = x # syntax sugar; macro-based logging @log input = u dx .= A*x + B*u end end
- (Optional)
Params(env::AbstractEnv)
: returns structured parameters of given environmentenv
.
Simulator
Simulator(state0, dyn, p; Problem, solver)
is a simulator struct that will be simulated bysolve
(non-interactive) orstep!
andstep_until!
(interactive).Problem = :ODE
andProblem = :Discrete
implyODEProblem
andDiscreteProblem
, respectively. For more details, seesrc/APIs/simulation.jl
.
Non-interactive interface (e.g., compatible with callbacks from DifferentialEquations.jl)
solve(simulation::Simulator)
will solve (O)DE and providedf::DataFrame
.- For now, only in-place method (iip) is supported.
Interactive interface (you should be aware of how to use integrator
interface in DifferentialEquations.jl)
reinit!(simulator::Simulator)
will reinitialisesimulator::Simulator
.step!(simulator::Simulator, Δt; stop_at_tdt=true)
will step thesimulator::Simulator
asΔt
.step_until!(simulator::Simulator, tf)
will step thesimulator::Simulator
untiltf
.push!(simulator::Simulator, df::DataFrame)
will push a datum fromsimulator
todf
.
Utilities
apply_inputs(func; kwargs...)
- By using this, user can easily apply external inputs into environments. It is borrowed from an MRAC example of ComponentArrays.jl and extended to be compatible with SimulationLogger.jl.
- (Limitations) for now, dynamical equations wrapped by
apply_inputs
will automatically generate logging function (even without@Loggable
). In this case, all data will be an array of emptyNamedTuple
.
- Macros for logging data:
@Loggable
,@log
,@onlylog
,@nested_log
,@nested_onlylog
.- For more details, see SimulationLogger.jl.
See directory ./test
.
using FSimBase
using DifferentialEquations
using ComponentArrays
using Test
using LinearAlgebra
using DataFrames
function main()
state0 = [1.0, 2.0]
p = 1
tf = 1.0
Δt = 0.01
@Loggable function dynamics!(dx, x, p, t)
@log t
@log x
dx .= -p.*x
end
simulator = Simulator(
state0, dynamics!, p;
Problem=ODEProblem,
solver=Tsit5(),
tf=tf,
)
# solve approach (automatically reinitialised)
@time _df = solve(simulator; savestep=Δt)
# interactive simulation
## step!
reinit!(simulator)
step!(simulator, Δt)
@test simulator.integrator.t ≈ Δt
## step_until! (callback-like)
ts_weird = 0:Δt:tf+Δt
df_ = DataFrame()
reinit!(simulator)
@time for t in ts_weird
flag = step_until!(simulator, t) # flag == false if step is inappropriate
if simulator.integrator.u[1] < 5e-1
break
else
push!(simulator, df_, flag) # push data only when flag == true
end
end
println(df_[end-5:end, :])
## step_until!
df = DataFrame()
reinit!(simulator)
@time for t in ts_weird
step_until!(simulator, t)
push!(simulator, df) # flag=true is default
# safer way:
# flag = step_until!(simulator, t)
# push!(simulator, df, flag)
# or, equivalently,
# push!(simulator, df, step_until!(simulator, t)) # compact form
end
println(df[end-5:end, :])
@test norm(_df.sol[end].x - df.sol[end].x) < 1e-6
@test simulator.integrator.t ≈ tf
end
@testset "minimal" begin
main()
end
- SimulationLogger.jl: A convenient logging tools compatible with DifferentialEquations.jl.