Steering

Description of the problem

The particle steering problem is a classical benchmark in time-optimal control. It models the task of steering a particle from a given initial state to a specified terminal state while minimising travel time. The system includes four state variables representing position and velocity components, and a single scalar control input that determines the steering direction. Control bounds and dynamics constraints must be satisfied throughout the trajectory [Athans & Falb 1966; Bryson & Ho 1975].

Mathematical formulation

The problem can be stated as

\[\begin{aligned} \min_{x,u} \quad & J(x,u) = t_f \\[0.5em] \text{s.t.} \quad & \dot{x}_1(t) = x_3(t), \\[0.5em] & \dot{x}_2(t) = x_4(t), \\[0.5em] & \dot{x}_3(t) = a \cos(u(t)), \\[0.5em] & \dot{x}_4(t) = a \sin(u(t)), \\[0.5em] & u_{\min} \le u(t) \le u_{\max}, \\[0.5em] & x(0) = x_s, \quad x(t_f) = x_f, \\[0.5em] & t_f \ge 0. \end{aligned}\]

where

  • the state $x = (x_1, x_2, x_3, x_4)$ are the particle states (position and velocity components),
  • the control $u$ is the steering control, bounded by $u_{\min} = -\pi/2$ and $u_{\max} = \pi/2$,
  • the parameter $a$ is a constant scaling the acceleration.

Parameter values

ParameterSymbolValue
Acceleration magnitude$a$100
Control lower bound$u_{\min}$$-\pi/2$
Control upper bound$u_{\max}$$\pi/2$
Initial state$x_s$$(0,0,0,0)$
Final state$x_f$$(\text{unspecified},5,45,0)$

Qualitative behaviour

The optimal trajectory typically exhibits bang–bang control, where the steering input switches between its extreme values to minimise time. The trajectory balances position and velocity to reach the terminal state in minimal time. Analytical solutions exist for simplified cases, while numerical methods are required for general conditions.

Characteristics

  • Four-dimensional state with single scalar control,
  • Time-optimal objective,
  • Control constraints and boundary conditions,
  • Bang–bang control structure with possible singular arcs,
  • Serves as a benchmark for time-optimal control methods.

References

  • Athans, M., & Falb, P.L. (1966). Optimal Control: An Introduction to the Theory and Its Applications. McGraw-Hill. This classic text introduces the theory of optimal control, including time-optimal problems and bang–bang solutions, providing foundational concepts relevant to particle steering problems.

  • Bryson, A.E., & Ho, Y.-C. (1975). Applied Optimal Control: Optimization, Estimation, and Control. Hemisphere Publishing. The book presents methods and examples of time-optimal and energy-optimal control problems, including systems with bang–bang and singular arc solutions, directly relevant to the steering problem.

  • More, J., Garbow, B., Hillstrom, K., & Watson, L. (2001). COPS: Constrained Optimization Problem Set (COPS3). Mathematics and Computer Science Division, Argonne National Laboratory. https://www.mcs.anl.gov/~more/cops/cops3.pdf This report lists the particle steering problem as a benchmark for nonlinear constrained optimization, providing the original formulation and serving as a reference for numerical experimentation with optimal control solvers.

  • PSOPT Example Set. Particle Steering Problem. https://www.psopt.net/list-of-examples?utm_source=chatgpt.com A modern implementation of the particle steering problem for testing direct optimal control methods, illustrating numerical approaches for time-optimal trajectory problems.

Numerical set-up

In this section, we prepare the numerical environment required to study the problem. We begin by importing the relevant Julia packages and then initialise the data frames that will store the results of our simulations and computations. These structures provide the foundation for solving the problem and for comparing the different solution strategies in a consistent way.

using OptimalControlProblems    # to access the Beam model
using OptimalControl            # to import the OptimalControl model
using NLPModelsIpopt            # to solve the model with Ipopt
import DataFrames: DataFrame    # to store data
using NLPModels                 # to retrieve data from the NLP solution
using Plots                     # to plot the trajectories
using Plots.PlotMeasures        # for leftmargin, bottommargin
using JuMP                      # to import the JuMP model
using Ipopt                     # to solve the JuMP model with Ipopt
using Printf                    # to print

data_pb = DataFrame(            # to store data about the problem
    Problem=Symbol[],
    Grid_Size=Int[],
    Variables=Int[],
    Constraints=Int[],
)

data_re = DataFrame(            # to store data about the resolutions
    Model=Symbol[],
    Flag=Any[],
    Iterations=Int[],
    Objective=Float64[],
)

Metadata

The default number of time steps is:

metadata(:steering)[:grid_size]
500

The default values of the parameters are:

metadata(:steering)[:parameters]
Parameter = Value
------------------
    t0 =  0.0000e+00
     a =  1.0000e+02
 u_min = -1.5708e+00
 u_max =  1.5708e+00
  tf_l =  0.0000e+00
 x₁_t0 =  0.0000e+00
 x₂_t0 =  0.0000e+00
 x₃_t0 =  0.0000e+00
 x₄_t0 =  0.0000e+00
 x₂_tf =  5.0000e+00
 x₃_tf =  4.5000e+01
 x₄_tf =  0.0000e+00

Initial guess

Before solving the problem, it is often useful to inspect the initial guess (sometimes called the first iterate). This guess is obtained by running the NLP solver with max_iter = 0, which evaluates the problem formulation without performing any optimisation steps.

We plot the resulting trajectories for both the OptimalControl and JuMP models. Since both backends represent the same mathematical problem, their initial guesses should coincide, providing a useful consistency check before moving on to the optimised solution.

Code to plot the initial guess
Click to unfold and see the code for plotting the initial guess.
function plot_initial_guess(problem)

    # -----------------------------
    # Build OptimalControl problem
    # -----------------------------
    docp = eval(problem)(OptimalControlBackend())
    nlp_oc = nlp_model(docp)
    ocp_oc = ocp_model(docp)

    # Solve NLP with zero iterations (initial guess)
    nlp_oc_sol = NLPModelsIpopt.ipopt(nlp_oc; max_iter=0)

    # Build OptimalControl solution
    ocp_sol = build_ocp_solution(docp, nlp_oc_sol)

    # get dimensions
    n = state_dimension(ocp_oc)
    m = control_dimension(ocp_oc)

    # -----------------------------
    # Plot OptimalControl solution
    # -----------------------------
    plt = plot(
        ocp_sol;
        state_style=(color=1,),
        costate_style=(color=1, legend=:none),
        control_style=(color=1, legend=:none),
        path_style=(color=1, legend=:none),
        dual_style=(color=1, legend=:none),
        size=(816, 220*(n+m)),
        label="OptimalControl",
        leftmargin=20mm,
    )

    # Hide legend for additional state plots
    for i in 2:n
        plot!(plt[i]; legend=:none)
    end

    # -----------------------------
    # Build JuMP model
    # -----------------------------
    nlp_jp = eval(problem)(JuMPBackend())

    # Solve NLP with zero iterations (initial guess)
    set_optimizer(nlp_jp, Ipopt.Optimizer)
    set_optimizer_attribute(nlp_jp, "max_iter", 0)
    optimize!(nlp_jp)

    # Extract trajectories
    t_grid = time_grid(nlp_jp)
    x_fun = state(nlp_jp)
    u_fun = control(nlp_jp)
    p_fun = costate(nlp_jp)

    # -----------------------------
    # Plot JuMP solution on top
    # -----------------------------
    # States
    for i in 1:n
        label = i == 1 ? "JuMP" : :none
        plot!(plt[i], t_grid, t -> x_fun(t)[i]; color=2, linestyle=:dash, label=label)
    end

    # Costates
    for i in 1:n
        plot!(plt[n+i], t_grid, t -> -p_fun(t)[i]; color=2, linestyle=:dash, label=:none)
    end

    # Controls
    for i in 1:m
        plot!(plt[2*n+i], t_grid, t -> u_fun(t)[i]; color=2, linestyle=:dash, label=:none)
    end

    return plt
end
plot_initial_guess(:steering)
Example block output

Solving the problem

To solve an optimal control problem, we can rely on two complementary formulations: the OptimalControl backend, which works directly with the discretised control problem, and the JuMP backend, which leverages JuMP’s flexible modelling framework.

Both approaches generate equivalent NLPs that can be solved with Ipopt, and comparing them ensures consistency between the two formulations.

Before solving, we can inspect the discretisation details of the problem. The table below reports the number of grid points, decision variables, and constraints associated with the chosen formulation.

push!(data_pb,(
    Problem=:steering,
    Grid_Size=metadata(:steering)[:grid_size],
    Variables=get_nvar(nlp_model(steering(OptimalControlBackend()))),
    Constraints=get_ncon(nlp_model(steering(OptimalControlBackend()))),
))
1×4 DataFrame
RowProblemGrid_SizeVariablesConstraints
SymbolInt64Int64Int64
1steering50025062007

OptimalControl model

We first solve the problem using the OptimalControl backend. The process begins by importing the problem definition and constructing the associated nonlinear programming (NLP) model. This NLP is then passed to the Ipopt solver, with standard options for tolerance and barrier parameter strategy.

# import DOCP model
docp = steering(OptimalControlBackend())

# get NLP model
nlp_oc = nlp_model(docp)

# solve
nlp_oc_sol = NLPModelsIpopt.ipopt(
    nlp_oc;
    print_level=4,
    tol=1e-8,
    mu_strategy="adaptive",
    sb="yes",
)
Total number of variables............................:     2506
                     variables with only lower bounds:        1
                variables with lower and upper bounds:      501
                     variables with only upper bounds:        0
Total number of equality constraints.................:     2007
Total number of inequality constraints...............:        0
        inequality constraints with only lower bounds:        0
   inequality constraints with lower and upper bounds:        0
        inequality constraints with only upper bounds:        0


Number of Iterations....: 17

                                   (scaled)                 (unscaled)
Objective...............:   5.5457186103692224e-01    5.5457186103692224e-01
Dual infeasibility......:   7.5784997368786306e-11    7.5784997368786306e-11
Constraint violation....:   3.5834613193452824e-09    3.5834613193452824e-09
Variable bound violation:   0.0000000000000000e+00    0.0000000000000000e+00
Complementarity.........:   1.0238142261644498e-11    1.0238142261644498e-11
Overall NLP error.......:   3.5834613193452824e-09    3.5834613193452824e-09


Number of objective function evaluations             = 19
Number of objective gradient evaluations             = 18
Number of equality constraint evaluations            = 19
Number of inequality constraint evaluations          = 0
Number of equality constraint Jacobian evaluations   = 18
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations             = 17
Total seconds in IPOPT                               = 0.633

EXIT: Optimal Solution Found.

JuMP model

We now repeat the procedure using the JuMP backend. Here, the problem is reformulated as a JuMP model, which offers a flexible and widely used framework for nonlinear optimisation in Julia. The solver settings are chosen to mirror those used previously, so that the results can be compared on an equal footing.

# import model
nlp_jp = steering(JuMPBackend())

# solve with Ipopt
set_optimizer(nlp_jp, Ipopt.Optimizer)
set_optimizer_attribute(nlp_jp, "print_level", 4)
set_optimizer_attribute(nlp_jp, "tol", 1e-8)
set_optimizer_attribute(nlp_jp, "mu_strategy", "adaptive")
set_optimizer_attribute(nlp_jp, "linear_solver", "mumps")
set_optimizer_attribute(nlp_jp, "sb", "yes")
optimize!(nlp_jp)
Total number of variables............................:     2506
                     variables with only lower bounds:        1
                variables with lower and upper bounds:      501
                     variables with only upper bounds:        0
Total number of equality constraints.................:     2007
Total number of inequality constraints...............:        0
        inequality constraints with only lower bounds:        0
   inequality constraints with lower and upper bounds:        0
        inequality constraints with only upper bounds:        0


Number of Iterations....: 17

                                   (scaled)                 (unscaled)
Objective...............:   5.5457186103692246e-01    5.5457186103692246e-01
Dual infeasibility......:   7.5782332833575182e-11    7.5782332833575182e-11
Constraint violation....:   3.5834614164897971e-09    3.5834614164897971e-09
Variable bound violation:   0.0000000000000000e+00    0.0000000000000000e+00
Complementarity.........:   1.0238142263371651e-11    1.0238142263371651e-11
Overall NLP error.......:   3.5834614164897971e-09    3.5834614164897971e-09


Number of objective function evaluations             = 19
Number of objective gradient evaluations             = 18
Number of equality constraint evaluations            = 19
Number of inequality constraint evaluations          = 0
Number of equality constraint Jacobian evaluations   = 18
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations             = 17
Total seconds in IPOPT                               = 0.136

EXIT: Optimal Solution Found.

Numerical comparisons

In this section, we examine the results of the problem resolutions. We extract the solver status (flag), the number of iterations, and the objective value for each model. This provides a first overview of how each approach performs and sets the stage for a more detailed comparison of the solution trajectories.

# from OptimalControl model
push!(data_re,(
    Model=:OptimalControl,
    Flag=nlp_oc_sol.status,
    Iterations=nlp_oc_sol.iter,
    Objective=nlp_oc_sol.objective,
))

# from JuMP model
push!(data_re,(
    Model=:JuMP,
    Flag=termination_status(nlp_jp),
    Iterations=barrier_iterations(nlp_jp),
    Objective=objective_value(nlp_jp),
))
2×4 DataFrame
RowModelFlagIterationsObjective
SymbolAnyInt64Float64
1OptimalControlfirst_order170.554572
2JuMPLOCALLY_SOLVED170.554572

We compare the solutions obtained from the OptimalControl and JuMP models by examining the number of iterations required for convergence, the $L^2$-norms of the differences in states, controls, and additional variables, and the corresponding objective values. Both absolute and relative errors are reported, providing a clear quantitative measure of the agreement between the two approaches.

Code to print the numerical comparisons
Click to unfold and get the code of the numerical comparisons.
function L2_norm(T, X)
    # T and X are supposed to be one dimensional
    s = 0.0
    for i in 1:(length(T) - 1)
        s += 0.5 * (X[i]^2 + X[i + 1]^2) * (T[i + 1]-T[i])
    end
    return √(s)
end

function print_numerical_comparisons(problem, docp, nlp_oc_sol, nlp_jp)

    # get relevant data from OptimalControl model
    ocp_sol = build_ocp_solution(docp, nlp_oc_sol)
    t_oc = time_grid(ocp_sol)
    x_oc = state(ocp_sol).(t_oc)
    u_oc = control(ocp_sol).(t_oc)
    v_oc = variable(ocp_sol)
    o_oc = objective(ocp_sol)
    i_oc = iterations(ocp_sol)

    # get relevant data from JuMP model
    t_jp = time_grid(nlp_jp)
    x_jp = state(nlp_jp).(t_jp)
    u_jp = control(nlp_jp).(t_jp)
    o_jp = objective(nlp_jp)
    v_jp = variable(nlp_jp)
    i_jp = iterations(nlp_jp)

    x_vars = state_components(nlp_jp)
    u_vars = control_components(nlp_jp)
    v_vars = variable_components(nlp_jp)

    println("┌─ ", string(problem))
    println("│")
    println("├─  Number of Iterations")
    @printf("│     OptimalControl : %d   JuMP : %d\n", i_oc, i_jp)

    # States
    println("├─  States (L2 Norms)")
    for i in eachindex(x_vars)
        xi_oc = [x_oc[k][i] for k in eachindex(t_oc)]
        xi_jp = [x_jp[k][i] for k in eachindex(t_jp)]
        L2_ae = L2_norm(t_oc, xi_oc - xi_jp)
        L2_re = L2_ae / (0.5 * (L2_norm(t_oc, xi_oc) + L2_norm(t_oc, xi_jp)))
        @printf("│     %-6s Abs: %.3e   Rel: %.3e\n", x_vars[i], L2_ae, L2_re)
    end

    # Controls
    println("├─  Controls (L2 Norms)")
    for i in eachindex(u_vars)
        ui_oc = [u_oc[k][i] for k in eachindex(t_oc)]
        ui_jp = [u_jp[k][i] for k in eachindex(t_jp)]
        L2_ae = L2_norm(t_oc, ui_oc - ui_jp)
        L2_re = L2_ae / (0.5 * (L2_norm(t_oc, ui_oc) + L2_norm(t_oc, ui_jp)))
        @printf("│     %-6s Abs: %.3e   Rel: %.3e\n", u_vars[i], L2_ae, L2_re)
    end

    # Variables
    if !isnothing(v_vars)
        println("├─  Variables")
        for i in eachindex(v_vars)
            vi_oc = v_oc[i]
            vi_jp = v_jp[i]
            vi_ae = abs(vi_oc - vi_jp)
            vi_re = vi_ae / (0.5 * (abs(vi_oc) + abs(vi_jp)))
            @printf("│     %-6s Abs: %.3e   Rel: %.3e\n", v_vars[i], vi_ae, vi_re)
        end
    end

    # Objective
    o_ae = abs(o_oc - o_jp)
    o_re = o_ae / (0.5 * (abs(o_oc) + abs(o_jp)))
    println("├─  Objective")
    @printf("│            Abs: %.3e   Rel: %.3e\n", o_ae, o_re)
    println("└─")
    return nothing
end
print_numerical_comparisons(:steering, docp, nlp_oc_sol, nlp_jp)
┌─ steering
│
├─  Number of Iterations
│     OptimalControl : 17   JuMP : 17
├─  States (L2 Norms)
│     x₁     Abs: 1.718e-15   Rel: 4.206e-16
│     x₂     Abs: 1.538e-15   Rel: 6.741e-16
│     x₃     Abs: 8.285e-15   Rel: 4.201e-16
│     x₄     Abs: 7.092e-15   Rel: 9.494e-16
├─  Controls (L2 Norms)
│     u      Abs: 7.419e-16   Rel: 1.581e-15
├─  Variables
│     tf     Abs: 2.220e-16   Rel: 4.004e-16
├─  Objective
│            Abs: 2.220e-16   Rel: 4.004e-16
└─

Plotting the solutions

In this section, we visualise the trajectories of the states, costates, and controls obtained from both the OptimalControl and JuMP solutions. The plots provide an intuitive way to compare the two approaches and to observe how the constraints and the optimal control influence the system dynamics.

For each variable, the OptimalControl solution is shown in solid lines, while the JuMP solution is overlaid using dashed lines. Since both models represent the same mathematical problem, their trajectories should closely coincide, highlighting the consistency between the two formulations.

# build an ocp solution to use the plot from OptimalControl package
ocp_sol = build_ocp_solution(docp, nlp_oc_sol)

# dimensions
n = state_dimension(ocp_sol)
m = control_dimension(ocp_sol)

# from OptimalControl solution
plt = plot(
    ocp_sol;
    color=1,
    size=(816, 240*(n+m)),
    label="OptimalControl",
    leftmargin=20mm,
)
for i in 2:length(plt)
    plot!(plt[i]; legend=:none)
end

# from JuMP solution
t = time_grid(nlp_jp)     # t0, ..., tN = tf
x = state(nlp_jp)         # function of time
u = control(nlp_jp)       # function of time
p = costate(nlp_jp)       # function of time

for i in 1:n # state
    label = i == 1 ? "JuMP" : :none
    plot!(plt[i], t, t -> x(t)[i]; color=2, linestyle=:dash, label=label)
end

for i in 1:n # costate
    plot!(plt[n+i], t, t -> -p(t)[i]; color=2, linestyle=:dash, label=:none)
end

for i in 1:m # control
    plot!(plt[2n+i], t, t -> u(t)[i]; color=2, linestyle=:dash, label=:none)
end
Example block output