FlightSims

FlightSims.jl is a general-purpose numerical simulator supporting nested environments and convenient macro-based data logging.

Road map

  • an easy routine for saving and loading data see DrWatson.jl for scientific projects.
  • ROS2 compatibility (not urgent)

Useful packages

Features

If you want more functionality, please feel free to report an issue!

Nested Environments and Zoo

  • Environments usually stand for dynamical systems but also include other utilities, for example, controllers.

  • One can generate user-defined nested environments using provided APIs. Also, some predefined environments are provided for reusability (i.e., environment zoo). Take a look at src/environments.

  • Examples include

    basics
    multicopters
    • (Hexacopter) LeeHexacopterEnv
    • (Quadcopter) IslamQuadcopterEnv, GoodarziQuadcopterEnv
    multicopters
    • (Hexacopter) LeeHexacopterEnv (currently maintained)
    • (Quadcopter) IslamQuadcopterEnv, GoodarziQuadcopterEnv
    allocators
    • (Moore-Penrose pseudo inverse control allocation) PseudoInverseAllocator
    controllers
    • (Linear quadratic regulator) LQR
    • (Proportional-Integral-Derivative controller) PID
      • Note that the derivative term is obtained via second-order filter.
    • (Pure proportional navigation guidance) PPNG
    integrated_environments
    • (Backstepping Position Controller + Static Allocator + Multicopter) BacksteppingPositionController_StaticAllocator_MulticopterEnv
      • For example, BacksteppingPositionControllerEnv (backstepping position controller) + PseudoInverseAllocator (pseudo-inverse allocator, a static allocator) + LeeHexacopterEnv (hexacopter, a multicopter)
    • See src/environments/integrated_environments.
    missiles
    • (point-mass simple missile in 3D space) PointMass3DMissile
    • (pursuer vs evador in 3D space) PursuerEvador3DMissile

Utilities

  • Some utilities are also provided for dynamical system simulation.
  • Examples include
    • Simulation rendering (currently maintained)
      • (Multicopter rendering) See src/environments/multicopters/render.jl.
    • Function approximator
      • (Approximator) LinearApproximator, PolynomialBasis
    • Data manipulation for machine learning
      • (Split data) partitionTrainTest
    • Reference trajectory generator
      • (Command generator) HelixCommandGenerator, PowerLoop

APIs

Main APIs are provided in src/APIs.

Make an environment

  • AbstractEnv: an abstract type for user-defined and predefined environments. In general, environments is a sub-type of AbstractEnv.
    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), Dynamics(env::AbstractEnv): return a function that maps in-place (recommended) and out-of-place dynamics (resp.), 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)
            @log state = x
            @log input = u
            dx .= A*x + B*u
        end
    end
    
  • (Optional) Params(env::AbstractEnv): returns parameters of given environment env.

Simulation, logging, and data saving & loading

Main APIs

  • sim
  • 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 empty NamedTuple.
  • Macros for logging data: @Loggable, @log, @onlylog, @nested_log, @nested_onlylog.
  • Example code
    A = [0 1;
         0 0]  # 2 x 2
    B = [0 1]'  # 2 x 1
    n, m = 2, 1
    env = LinearSystemEnv(A, B)  # exported from FlightSims
    x0 = State(env)([1.0, 2.0])
    # optimal control
    Q = Matrix(I, n, n)
    R = Matrix(I, m, m)
    lqr = LQR(A, B, Q, R)  # exported from FlightSims
    u_lqr = FS.OptimalController(lqr)  # (x, p, t) -> -K*x; minimise J = ∫ (x' Q x + u' R u) from 0 to ∞
    
    # simulation
    tf = 10.0
    Δt = 0.01
    prob, df = sim(
                   x0,  # initial condition
                   apply_inputs(dynamics!; u=(x, p, t) -> u_lqr(x));
                   tf=tf,  # final time
                   savestep=Δt,
                  )
    

Examples

Basic: toy example

  • For an introductory example of FlightSims.jl, see the following example code (test/toy_example.jl).
using FlightSims
const FS = FlightSims

using LinearAlgebra  # for I, e.g., Matrix(I, n, n)
using ComponentArrays
using UnPack
using Transducers
using Plots
using DifferentialEquations  # for callbacks


struct MyEnv <: AbstractEnv  # AbstractEnv exported from FS
    a
    b
end

"""
FlightSims recommends you to use closures for State and Dynamics!. For more details, see https://docs.julialang.org/en/v1/devdocs/functions/.
"""
function State(env::MyEnv)
    return function (x1::Number, x2::Number)
        ComponentArray(x1=x1, x2=x2)
    end
end

function Dynamics!(env::MyEnv)
    @unpack a, b = env  # @unpack is very useful!
    @Loggable function dynamics!(dx, x, p, t; u)  # `Loggable` makes it loggable via SimulationLogger.jl (imported in FS)
        @unpack x1, x2 = x
        @log x1  # to log x1
        @log x2  # to log x2
        dx.x1 = a*x2
        dx.x2 = b*u
    end
end

function my_controller(x, p, t)
    @unpack x1, x2 = x
    -(x1+x2)
end


function main()
    n = 2
    m = 1
    a, b = 1, 1
    A = -Matrix(I, n, n)
    B = Matrix(I, m, m)
    env = MyEnv(a, b)
    tf = 10.0
    Δt = 0.01
    # callbacks
    function condition(u, t, integrator)  # DiffEq.jl convention
        @unpack x1, x2 = u
        x1 - 0  # i.e., it stops when x1 = 0
    end
    affect!(integrator) = terminate!(integrator)  # See DiffEq.jl documentation
    cb_diverge = ContinuousCallback(condition, affect!)
    cb = CallbackSet(cb_diverge,)  # useful for multiple callbacks
    x10, x20 = 10.0, 0.0
    x0 = State(env)(x10, x20)
    # prob: DE problem, df: DataFrame
    @time prob, df = sim(
                         x0,  # initial condition
                         apply_inputs(Dynamics!(env); u=my_controller);  # dynamics!; apply_inputs is exported from FS and is so useful for systems with inputs
                         tf=10.0,
                         savestep=Δt,  # savestep is NOT simulation step
                         callback=cb,
                        )  # sim is exported from FS
    ts = df.time
    x1s = df.sol |> Map(datum -> datum.x1) |> collect
    x2s = df.sol |> Map(datum -> datum.x2) |> collect
    # plot
    p_x1 = plot(ts, x1s;
                label="x1",
               )
    p_x2 = plot(ts, x2s;
                label="x2",
               )
    p_x = plot(p_x1, p_x2, layout=(2, 1))
    # save
    dir_log = "figures"
    mkpath(dir_log)
    savefig(p_x, joinpath(dir_log, "toy_example.png"))
    display(p_x)
end

ex_screenshot

Optimal control and reinforcement learning

  • For an example of infinite-horizon continuous-time linear quadratic regulator (LQR), see the following example code (test/lqr.jl).
using FlightSims
const FS = FlightSims
using DifferentialEquations
using LinearAlgebra
using Plots
using Test
using Transducers


function test()
    # linear system
    A = [0 1;
         0 0]  # 2 x 2
    B = [0 1]'  # 2 x 1
    n, m = 2, 1
    env = LinearSystemEnv(A, B)  # exported from FlightSims
    x0 = State(env)([1.0, 2.0])
    p0 = zero.(x0)  # auxiliary parameter
    # optimal control
    Q = Matrix(I, n, n)
    R = Matrix(I, m, m)
    lqr = LQR(A, B, Q, R)  # exported from FlightSims
    u_lqr = FS.OptimalController(lqr)  # (x, p, t) -> -K*x; minimise J = ∫ (x' Q x + u' R u) from 0 to ∞

    # simulation
    tf = 10.0
    Δt = 0.01
    affect!(integrator) = integrator.p = copy(integrator.u)  # auxiliary callback funciton
    cb = PeriodicCallback(affect!, Δt; initial_affect=true)  # auxiliary callback
    @Loggable function dynamics!(dx, x, p, t)
        @onlylog p  # activate this line only when logging data
        u = u_lqr(x)
        @log x, u
        @nested_log Dynamics!(env)(dx, x, p, t; u=u)  # exported `state` and `input` from `Dynamics!(env)`
    end
    prob, df = sim(
                   x0,  # initial condition
                   dynamics!,  # dynamics with input of LQR
                   p0;
                   tf=tf,  # final time
                   callback=cb,
                   savestep=Δt,
                  )
    ts = df.time
    xs = df.sol |> Map(datum -> datum.x) |> collect
    us = df.sol |> Map(datum -> datum.u) |> collect
    ps = df.sol |> Map(datum -> datum.p) |> collect
    states = df.sol |> Map(datum -> datum.state) |> collect
    inputs = df.sol |> Map(datum -> datum.input) |> collect
    @test xs == states
    @test us == inputs
    p_x = plot(ts, hcat(states...)';
               title="state variable", label=["x1" "x2"], color=[:black :black], lw=1.5,
              )  # Plots
    plot!(p_x, ts, hcat(ps...)';
          ls=:dash, label="param", color=[:red :orange], lw=1.5
         )
    savefig("figures/x_lqr.png")
    plot(ts, hcat(inputs...)'; title="control input", label="u")  # Plots
    savefig("figures/u_lqr.png")
    df
end
julia> test()
1001×2 DataFrame
  Row  time     sol
       Float64  NamedTup
──────┼────────────────────────────────────────────
    1     0.0   (p = [1.01978, 1.95564], state =…
    2     0.01  (p = [1.01978, 1.95564], state =…
    3     0.02  (p = [1.03911, 1.91186], state =…
    4     0.03  (p = [1.05802, 1.86863], state =…
    5     0.04  (p = [1.07649, 1.82596], state =…
                              
  998     9.97  (p = [-0.00093419, 0.00103198], 
  999     9.98  (p = [-0.000923913, 0.00102347],
 1000     9.99  (p = [-0.00091372, 0.001015], st
 1001    10.0   (p = [-0.00091372, 0.001015], st
                                   992 rows omitted

ex_screenshot ex_screenshot

Multicopter position control

  • For an example of backstepping position tracking controller for quadcopters, see test/environments/integrated_environments/backstepping_position_controller_static_allocator_multicopter_env.jl.

Missile guidance with interactive visualisation

  • See test/pluto_guidance.jl (thanks to @nhcho91).

Alt Text

Multicopter rendering

  • See test/render.jl.

Alt Text

Visualisation of hexacopter including reference frame and topview

ex_screenshot ex_screenshot