Rocket

Description of the problem

The Goddard rocket problem is a classical optimal control problem that models the ascent of a vertically launched rocket. The objective is to maximise the final altitude by optimally controlling the thrust, while accounting for atmospheric drag, gravity, and fuel consumption. This problem is notable for its inclusion of singular arcs, where the control is neither at its maximum nor minimum, complicating the solution process. Historically, it was first proposed by Robert H. Goddard in 1919 (A Method of Reaching Extreme Altitudes), one of the founding works in modern rocketry.

Mathematical formulation

The problem can be stated as

\[\begin{aligned} \min_{h,v,m,T} \quad & J = -h(T) \\[0.5em] \text{s.t.} \quad \dot{h}(t) &= v(t), \\[0.25em] \dot{v}(t) &= \frac{T(t) - D(h(t), v(t)) - m(t) g(h(t))}{m(t)}, \\[0.25em] \dot{m}(t) &= -\frac{T(t)}{c}, \\[0.25em] 0 \le T(t) &\le T_{\max}, \\[0.25em] h(t) &\ge h_0, \\[0.25em] v(t) &\ge v_0, \\[0.25em] m_f \le m(t) &\le m_0, \\[0.25em] h(0) &= h_0, \quad v(0) = v_0, \quad m(0) = m_0, \\[0.25em] m(T) &= m_f, \end{aligned}\]

where

  • Drag:

\[D(h,v) = D_c \, v^2 \exp\!\left(-h_c \frac{h - h_0}{h_0}\right)\]

  • Gravity:

\[g(h) = g_0 \left(\frac{h_0}{h}\right)^2\]

  • Fuel constant:

\[c = \frac{1}{2} \sqrt{g_0 h_0}\]

  • Maximum thrust:

\[T_{\max} = T_c m_0 g_0\]

Parameters

ParameterSymbolValue
Initial altitude$h_0$1
Initial velocity$v_0$0
Initial mass$m_0$1
Gravitational constant$g_0$1
Thrust coefficient$T_c$3.5
Drag coefficient$D_c$$\frac{1}{2} v_c \frac{m_0}{g_0}$
Characteristic altitude$h_c$500
Characteristic velocity$v_c$620
Characteristic mass ratio$m_c$0.6
Final mass$m_f$$m_c m_0$

Qualitative behaviour

The optimal trajectory often exhibits a bang–singular–bang structure, with thrust at maximum, then along a singular arc, and finally at minimum. The singular arc arises from the interplay of drag and gravity, producing a non-trivial optimal control law. This behaviour illustrates the complexity of the problem and the challenges in numerically resolving singular arcs.

Characteristics

  • Nonlinear coupled dynamics with drag and gravity,
  • Singular arcs in the optimal control,
  • State and control constraints,
  • Classical benchmark for testing optimal control algorithms, especially in aerospace applications.

References

  • Goddard, R. H. (1919). A Method of Reaching Extreme Altitudes. Smithsonian Institution, Washington D.C. The original formulation of the vertical rocket ascent problem, introducing fuel consumption, gravity, and drag effects.

  • Bryson, A. E. (1999). Dynamic Optimization. Addison Wesley Longman. pp. 392–394. Provides an introduction to singular arcs in optimal control with the Goddard rocket as a primary example.

  • Garfinkel, B. (1963). A Solution of the Goddard Problem. SIAM Journal on Control, 1(1), 1–20. doi:10.1137/0301020 Analytical and numerical solutions of the Goddard problem highlighting singular arc structures.

  • Tsiotras, P., & Kelley, H. J. (1992). Goddard Problem with Constrained Time of Flight. Journal of Guidance, Control, and Dynamics, 15(2), 394–399. doi:10.2514/3.20836 Study of time-constrained versions of the Goddard problem and associated optimal control strategies.

  • Bonnans, F., Martinon, P., & Trélat, E. (2007). Singular Arcs in the Generalized Goddard's Problem. arXiv:math/0703911 Analysis of singular arcs and numerical methods for the Goddard problem.

  • More, J., Garbow, B., Hillstrom, K., & Watson, L. (2001). COPS: Constrained Optimization Problem Set (COPS3). Argonne National Laboratory. Link Provides a benchmark implementation of the Goddard rocket problem for testing numerical optimization methods.

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(:rocket)[:grid_size]
500

The default values of the parameters are:

metadata(:rocket)[:parameters]
Parameter = Value
------------------
    t0 =  0.0000e+00
  h_t0 =  1.0000e+00
  v_t0 =  0.0000e+00
  m_t0 =  1.0000e+00
    g0 =  1.0000e+00
    Tc =  3.5000e+00
    hc =  5.0000e+02
    vc =  6.2000e+02
    mc =  6.0000e-01
   T_l =  0.0000e+00
  tf_l =  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(:rocket)
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=:rocket,
    Grid_Size=metadata(:rocket)[:grid_size],
    Variables=get_nvar(nlp_model(rocket(OptimalControlBackend()))),
    Constraints=get_ncon(nlp_model(rocket(OptimalControlBackend()))),
))
1×4 DataFrame
RowProblemGrid_SizeVariablesConstraints
SymbolInt64Int64Int64
1rocket50020051504

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 = rocket(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............................:     2005
                     variables with only lower bounds:     1003
                variables with lower and upper bounds:     1002
                     variables with only upper bounds:        0
Total number of equality constraints.................:     1504
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....: 21

                                   (scaled)                 (unscaled)
Objective...............:  -1.0128367166830137e+00    1.0128367166830137e+00
Dual infeasibility......:   1.3598075163735353e-09    1.3598075163735353e-09
Constraint violation....:   1.5371388051299562e-09    1.5371388051299562e-09
Variable bound violation:   1.1135572325331747e-37    1.1135572325331747e-37
Complementarity.........:   1.9466553661323338e-11    1.9466553661323338e-11
Overall NLP error.......:   1.5371388051299562e-09    1.5371388051299562e-09


Number of objective function evaluations             = 22
Number of objective gradient evaluations             = 22
Number of equality constraint evaluations            = 22
Number of inequality constraint evaluations          = 0
Number of equality constraint Jacobian evaluations   = 22
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations             = 21
Total seconds in IPOPT                               = 0.826

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 = rocket(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............................:     2005
                     variables with only lower bounds:     1003
                variables with lower and upper bounds:     1002
                     variables with only upper bounds:        0
Total number of equality constraints.................:     1504
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....: 21

                                   (scaled)                 (unscaled)
Objective...............:  -1.0128367166830152e+00    1.0128367166830152e+00
Dual infeasibility......:   1.3598087584345259e-09    1.3598087584345259e-09
Constraint violation....:   1.5371387773743805e-09    1.5371387773743805e-09
Variable bound violation:   4.6478574930930065e-36    4.6478574930930065e-36
Complementarity.........:   1.9466553740652001e-11    1.9466553740652001e-11
Overall NLP error.......:   1.5371387773743805e-09    1.5371387773743805e-09


Number of objective function evaluations             = 22
Number of objective gradient evaluations             = 22
Number of equality constraint evaluations            = 22
Number of inequality constraint evaluations          = 0
Number of equality constraint Jacobian evaluations   = 22
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations             = 21
Total seconds in IPOPT                               = 0.233

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_order211.01284
2JuMPLOCALLY_SOLVED211.01284

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(:rocket, docp, nlp_oc_sol, nlp_jp)
┌─ rocket
│
├─  Number of Iterations
│     OptimalControl : 21   JuMP : 21
├─  States (L2 Norms)
│     h      Abs: 6.021e-16   Rel: 1.341e-15
│     v      Abs: 8.431e-15   Rel: 2.564e-13
│     m      Abs: 1.328e-14   Rel: 4.396e-14
├─  Controls (L2 Norms)
│     T      Abs: 9.003e-12   Rel: 1.188e-11
├─  Variables
│     tf     Abs: 8.354e-15   Rel: 4.201e-14
├─  Objective
│            Abs: 1.554e-15   Rel: 1.535e-15
└─

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