Moonlander

This problem models the optimal control of a lunar lander aiming to reach a target position on the lunar surface while minimizing fuel consumption. The system is described by a set of differential equations governing the motion of the lander, including its position, velocity, orientation, and angular velocity.

Problem Formulation

The objective is to minimize the final time $t_f$ subject to the following constraints and dynamics:

  • State Variables: Position ($p_1$, $p_2$), velocity ($dp_1$, $dp_2$), orientation ($\theta$), and angular velocity ($d\theta$).
  • Control Inputs: Thrust forces $F_1$ and $F_2$ applied along the lander's orientation.
  • Final Time Constraint: $0.1 \le t_f \le 1.0$.
  • Control Constraints: $0 \le F_1(t), F_2(t) \le \text{max\_thrust}$.
  • Initial Conditions:

\[p_1(0) = 0, \quad p_2(0) = 0, \quad dp_1(0) = 0, \quad dp_2(0) = 0, \quad \theta(0) = 0, \quad d\theta(0) = 0\]

  • Final Conditions:

\[p_1(t_f) = \text{target}[1], \quad p_2(t_f) = \text{target}[2], \quad dp_1(t_f) = 0, \quad dp_2(t_f) = 0\]

System Dynamics

The dynamics of the lander are governed by

\[\dot{p}_1 = dp_1\]

\[\dot{p}_2 = dp_2\]

\[\dot{dp}_1 = \frac{1}{m} F_{\rm tot}[1]\]

\[\dot{dp}_2 = \frac{1}{m} F_{\rm tot}[2] - g\]

\[\dot{\theta} = d\theta\]

\[\dot{d\theta} = \frac{1}{I} \left( \frac{D}{2} (F_2 - F_1) \right)\]

where

\[F_{\rm tot} = \begin{bmatrix} \cos\theta & -\sin\theta \\ \sin\theta & \cos\theta \end{bmatrix} \begin{bmatrix} 0 \\ F_1 + F_2 \end{bmatrix}_{1:2}\]

represents the total thrust vector in the inertial frame, $m$ is the mass of the lander, $I$ its moment of inertia, $D$ the distance between thrusters, and $g$ the lunar gravitational acceleration.

Objective

The goal is to minimize the final time $t_f$:

\[J = t_f \to \min\]

subject to the dynamics, boundary conditions, and control constraints.

References

  • Vanroye, L., Sathya, A., De Schutter, J., & Decré, W. (2023). FATROP: A Fast Constrained Optimal Control Problem Solver for Robot Trajectory Optimization and Control. arXiv preprint arXiv:2303.16746. Retrieved from https://arxiv.org/pdf/2303.16746
  • Bryson, A. E., & Ho, Y.-C. (1975). Applied Optimal Control: Optimization, Estimation, and Control. Hemisphere Publishing.
  • Hull, D. G. (2007). Fundamentals of Airplane Flight Mechanics (Chapter on spacecraft and landing). Springer.

Packages

Import all necessary packages and define DataFrames to store information about the problem and resolution results.

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

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[],
)

Initial guess

The initial guess (or first iterate) can be visualised by running the solver with max_iter=0. Here is the initial guess.

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

    # dimensions
    x_vars = metadata[problem][:state_name]
    u_vars = metadata[problem][:control_name]
    n = length(x_vars) # number of states
    m = length(u_vars) # number of controls

    # import OptimalControl model
    docp = eval(problem)(OptimalControlBackend())
    nlp_oc = nlp_model(docp)

    # solve
    nlp_oc_sol = NLPModelsIpopt.ipopt(nlp_oc; max_iter=0)

    # build an optimal control solution
    ocp_sol = build_ocp_solution(docp, nlp_oc_sol)

    # plot the 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,
    )
    for i in 2:n
        plot!(plt[i]; legend=:none)
    end

    # import JuMP model
    nlp_jp = eval(problem)(JuMPBackend())

    # solve
    set_optimizer(nlp_jp, Ipopt.Optimizer)
    set_optimizer_attribute(nlp_jp, "max_iter", 0)
    optimize!(nlp_jp)

    # plot
    t = time_grid(problem, nlp_jp)     # t0, ..., tN = tf
    x = state(problem, nlp_jp)         # function of time
    u = control(problem, nlp_jp)       # function of time
    p = costate(problem, 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

    return plt
end

plot_initial_guess(:moonlander)
Example block output

Solve the problem

OptimalControl model

Import the OptimalControl model and solve it.

# import DOCP model
docp = moonlander(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............................:     4009
                     variables with only lower bounds:        0
                variables with lower and upper bounds:     1003
                     variables with only upper bounds:        0
Total number of equality constraints.................:     3010
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....: 52

                                   (scaled)                 (unscaled)
Objective...............:   9.6181005310230261e-01    9.6181005310230261e-01
Dual infeasibility......:   2.6617528320338479e-10    2.6617528320338479e-10
Constraint violation....:   2.9204114682102045e-09    2.9204114682102045e-09
Variable bound violation:   1.0855984555746545e-07    1.0855984555746545e-07
Complementarity.........:   1.6476450556781885e-11    1.6476450556781885e-11
Overall NLP error.......:   2.9204114682102045e-09    2.9204114682102045e-09


Number of objective function evaluations             = 72
Number of objective gradient evaluations             = 53
Number of equality constraint evaluations            = 72
Number of inequality constraint evaluations          = 0
Number of equality constraint Jacobian evaluations   = 53
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations             = 52
Total seconds in IPOPT                               = 3.066

EXIT: Optimal Solution Found.

The problem has the following numbers of steps, variables and constraints.

push!(data_pb,
    (
        Problem=:moonlander,
        Grid_Size=metadata[:moonlander][:N],
        Variables=get_nvar(nlp_oc),
        Constraints=get_ncon(nlp_oc),
    )
)
1×4 DataFrame
RowProblemGrid_SizeVariablesConstraints
SymbolInt64Int64Int64
1moonlander50040093010

JuMP model

Import the JuMP model and solve it.

# import model
nlp_jp = moonlander(JuMPBackend())

# solve
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............................:     4009
                     variables with only lower bounds:        0
                variables with lower and upper bounds:     1003
                     variables with only upper bounds:        0
Total number of equality constraints.................:     3010
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....: 46

                                   (scaled)                 (unscaled)
Objective...............:   9.6181005310879608e-01    9.6181005310879608e-01
Dual infeasibility......:   4.4554157141007272e-10    4.4554157141007272e-10
Constraint violation....:   6.0742708640759702e-09    6.0742708640759702e-09
Variable bound violation:   1.0855999832415364e-07    1.0855999832415364e-07
Complementarity.........:   2.3741629783139485e-11    2.3741629783139485e-11
Overall NLP error.......:   6.0742708640759702e-09    6.0742708640759702e-09


Number of objective function evaluations             = 70
Number of objective gradient evaluations             = 47
Number of equality constraint evaluations            = 70
Number of inequality constraint evaluations          = 0
Number of equality constraint Jacobian evaluations   = 47
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations             = 46
Total seconds in IPOPT                               = 0.956

EXIT: Optimal Solution Found.

Numerical comparisons

Let's get the flag, the number of iterations and the objective value from the resolutions.

# 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_order520.96181
2JuMPLOCALLY_SOLVED460.96181

We compare the OptimalControl and JuMP solutions in terms of the number of iterations, the $L^2$-norm of the differences in the state, control, and variable, as well as the objective values. Both absolute and relative errors are reported.

Click to unfold and get the code of the numerical comparison.
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 numerical_comparison(problem, docp, nlp_oc_sol, nlp_jp)

    # get relevant data from OptimalControl model
    ocp_sol = build_ocp_solution(docp, nlp_oc_sol) # build an ocp solution
    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(problem, nlp_jp)
    x_jp = state(problem, nlp_jp).(t_jp)
    u_jp = control(problem, nlp_jp).(t_jp)
    o_jp = objective(problem, nlp_jp)
    v_jp = variable(problem, nlp_jp)
    i_jp = iterations(problem, nlp_jp)

    x_vars = metadata[problem][:state_name]
    u_vars = metadata[problem][:control_name]
    v_vars = metadata[problem][:variable_name]

    println("┌─ ", string(problem))
    println("│")

    # number of iterations
    println("├─  Number of iterations")
    println("│")
    println("│     OptimalControl : ", i_oc)
    println("│     JuMP           : ", i_jp)
    println("│")

    # state
    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_oc = L2_norm(t_oc, xi_oc)
        L2_jp = L2_norm(t_oc, xi_jp)
        L2_ae = L2_norm(t_oc, xi_oc-xi_jp)
        L2_re = L2_ae/(0.5*(L2_oc + L2_jp))

        println("├─  State $(x_vars[i]) (L2 norm)")
        println("│")
        #println("│     OptimalControl : ", L2_oc)
        #println("│     JuMP           : ", L2_jp)
        println("│     Absolute error : ", L2_ae)
        println("│     Relative error : ", L2_re)
        println("│")
    end

    # control
    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_oc = L2_norm(t_oc, ui_oc)
        L2_jp = L2_norm(t_oc, ui_jp)
        L2_ae = L2_norm(t_oc, ui_oc-ui_jp)
        L2_re = L2_ae/(0.5*(L2_oc + L2_jp))

        println("├─  Control $(u_vars[i]) (L2 norm)")
        println("│")
        #println("│     OptimalControl : ", L2_oc)
        #println("│     JuMP           : ", L2_jp)
        println("│     Absolute error : ", L2_ae)
        println("│     Relative error : ", L2_re)
        println("│")
    end

    # variable
    if !isnothing(v_vars)
        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)))

            println("├─  Variable $(v_vars[i])")
            println("│")
            #println("│     OptimalControl : ", vi_oc)
            #println("│     JuMP           : ", vi_jp)
            println("│     Absolute error : ", vi_ae)
            println("│     Relative error : ", vi_re)
            println("│")
        end
    end

    # objective
    o_ae = abs(o_oc-o_jp)
    o_re = o_ae/(0.5*(abs(o_oc) + abs(o_jp)))

    println("├─  objective")
    println("│")
    #println("│     OptimalControl : ", o_oc)
    #println("│     JuMP           : ", o_jp)
    println("│     Absolute error : ", o_ae)
    println("│     Relative error : ", o_re)
    println("│")
    println("└─")

    return nothing
end

numerical_comparison(:moonlander, docp, nlp_oc_sol, nlp_jp)
┌─ moonlander
│
├─  Number of iterations
│
│     OptimalControl : 52
│     JuMP           : 46
│
├─  State p1 (L2 norm)
│
│     Absolute error : 1.3236938720326307e-7
│     Relative error : 4.702066706338599e-8
│
├─  State p2 (L2 norm)
│
│     Absolute error : 3.682687066729244e-7
│     Relative error : 1.2862853254481652e-7
│
├─  State dp1 (L2 norm)
│
│     Absolute error : 9.75598131544767e-6
│     Relative error : 1.5780866418218059e-6
│
├─  State dp2 (L2 norm)
│
│     Absolute error : 7.311375142926814e-6
│     Relative error : 1.3054722094536894e-6
│
├─  State θ (L2 norm)
│
│     Absolute error : 8.716662876024948e-7
│     Relative error : 6.857312706779865e-7
│
├─  State dθ (L2 norm)
│
│     Absolute error : 6.028772795420415e-5
│     Relative error : 9.708221576349734e-6
│
├─  Control F1 (L2 norm)
│
│     Absolute error : 0.11431887704631577
│     Relative error : 0.006788467801673329
│
├─  Control F2 (L2 norm)
│
│     Absolute error : 3.0177778937638107e-5
│     Relative error : 1.7458181059813582e-6
│
├─  Variable tf
│
│     Absolute error : 6.4934724264276156e-12
│     Relative error : 6.751304382254175e-12
│
├─  objective
│
│     Absolute error : 6.4934724264276156e-12
│     Relative error : 6.751304382254175e-12
│
└─

Plot the solutions

Visualise states, costates, and controls from the OptimalControl and JuMP solutions:

# 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)   # or length(metadata[:moonlander][:state_name])
m = control_dimension(ocp_sol) # or length(metadata[:moonlander][:control_name])

# from 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, 240*(n+m)),
    label="OptimalControl",
    leftmargin=20mm,
)
for i in 2:n
    plot!(plt[i]; legend=:none)
end

# from JuMP solution
t = time_grid(:moonlander, nlp_jp)     # t0, ..., tN = tf
x = state(:moonlander, nlp_jp)         # function of time
u = control(:moonlander, nlp_jp)       # function of time
p = costate(:moonlander, 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