Direct and indirect methods for the Goddard problem

Introduction

For this example, we consider the well-known Goddard problem[1] [2] which models the ascent of a rocket through the atmosphere, and we restrict here ourselves to vertical (one dimensional) trajectories. The state variables are the altitude $r$, speed $v$ and mass $m$ of the rocket during the flight, for a total dimension of 3. The rocket is subject to gravity $g$, thrust $u$ and drag force $D$ (function of speed and altitude). The final time $t_f$ is free, and the objective is to reach a maximal altitude with a bounded fuel consumption.

We thus want to solve the optimal control problem in Mayer form

\[ r(t_f) \to \max\]

subject to the controlled dynamics

\[ \dot{r} = v, \quad \dot{v} = \frac{T_{\max}\,u - D(r,v)}{m} - g, \quad \dot{m} = -u,\]

and subject to the control constraint $u(t) \in [0,1]$ and the state constraint $v(t) \leq v_{\max}$. The initial state is fixed while only the final mass is prescribed.

Nota bene

The Hamiltonian is affine with respect to the control, so singular arcs may occur, as well as constrained arcs due to the path constraint on the velocity (see below).

We import the OptimalControl.jl package to define the optimal control problem and NLPModelsIpopt.jl to solve it. We import the Plots.jl package to plot the solution. The OrdinaryDiffEq.jl package is used to define the shooting function for the indirect method and the NonlinearSolve.jl and MINPACK.jl packages permit to solve the shooting equation.

using OptimalControl  # to define the optimal control problem and more
using NLPModelsIpopt  # to solve the problem via a direct method
using OrdinaryDiffEq  # to get the Flow function from OptimalControl
using NonlinearSolve  # interface to NLE solvers
using MINPACK         # NLE solver: use to solve the shooting equation
using Plots           # to plot the solution

Optimal control problem

We define the problem

t0 = 0      # initial time
r0 = 1      # initial altitude
v0 = 0      # initial speed
m0 = 1      # initial mass
vmax = 0.1  # maximal authorized speed
mf = 0.6    # final mass to target

ocp = @def begin # definition of the optimal control problem

    tf ∈ R, variable
    t ∈ [t0, tf], time
    x = (r, v, m) ∈ R³, state
    u ∈ R, control

    x(t0) == [ r0, v0, m0 ]
    m(tf) == mf,         (1)
    0 ≤ u(t) ≤ 1
    r(t) ≥ r0
    0 ≤ v(t) ≤ vmax

    ẋ(t) == F0(x(t)) + u(t) * F1(x(t))

    r(tf) → max

end;

# Dynamics
const Cd = 310
const Tmax = 3.5
const β = 500
const b = 2

F0(x) = begin
    r, v, m = x
    D = Cd * v^2 * exp(-β*(r - 1)) # Drag force
    return [ v, -D/m - 1/r^2, 0 ]
end

F1(x) = begin
    r, v, m = x
    return [ 0, Tmax/m, -b*Tmax ]
end

Direct method

We then solve it

direct_sol = solve(ocp; grid_size=100)
This is Ipopt version 3.14.17, running with linear solver MUMPS 5.7.3.

Number of nonzeros in equality constraint Jacobian...:     2104
Number of nonzeros in inequality constraint Jacobian.:        0
Number of nonzeros in Lagrangian Hessian.............:     1111

Total number of variables............................:      405
                     variables with only lower bounds:      101
                variables with lower and upper bounds:      202
                     variables with only upper bounds:        0
Total number of equality constraints.................:      304
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

iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
   0 -1.0100000e+00 9.00e-01 2.00e+00   0.0 0.00e+00    -  0.00e+00 0.00e+00   0
   1 -1.0090670e+00 8.99e-01 6.67e+01   1.3 1.67e+02    -  3.64e-03 5.93e-04f  1
   2 -1.0000907e+00 8.74e-01 1.83e+02   1.0 6.64e+00    -  3.63e-02 2.83e-02h  1
   3 -1.0023670e+00 8.37e-01 1.34e+04   1.0 6.91e+00    -  2.11e-01 4.19e-02f  1
   4 -1.0025025e+00 7.70e-01 9.45e+03   1.5 4.04e+00    -  1.00e+00 8.09e-02f  1
   5 -1.0033626e+00 7.16e-01 1.48e+05   2.3 3.56e+00    -  3.73e-01 6.94e-02f  1
   6 -1.0142503e+00 9.62e-03 3.99e+04   2.3 7.16e-01    -  4.49e-01 9.90e-01h  1
   7 -1.0101264e+00 4.21e-03 4.24e+05   1.8 5.32e-01    -  5.03e-01 9.90e-01h  1
   8 -1.0068427e+00 3.20e-04 2.87e+06   0.9 2.44e-01    -  6.73e-01 9.91e-01h  1
   9 -1.0067336e+00 2.64e-06 2.30e+07   0.1 7.39e-02    -  7.07e-01 1.00e+00f  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  10 -1.0067340e+00 1.13e-10 6.50e+05  -5.0 2.90e-04    -  9.89e-01 1.00e+00h  1
  11 -1.0067350e+00 2.81e-10 7.20e+03  -7.0 4.26e-04    -  9.89e-01 1.00e+00h  1
  12 -1.0078967e+00 9.07e-04 5.95e+03  -3.0 7.61e-01    -  6.55e-01 7.21e-01f  1
  13 -1.0081866e+00 3.67e-06 9.38e+03  -9.0 1.31e-02    -  8.60e-01 1.00e+00h  1
  14 -1.0091814e+00 1.57e-04 1.79e+02  -4.6 2.18e-01    -  1.00e+00 9.28e-01h  1
  15 -1.0105115e+00 2.55e-04 7.17e+02  -4.4 3.09e-01    -  1.00e+00 5.75e-01h  1
  16 -1.0114149e+00 2.34e-05 7.98e-04  -5.1 4.29e-02    -  1.00e+00 1.00e+00h  1
  17 -1.0122891e+00 7.90e-05 1.24e+02  -5.8 1.25e-01    -  9.98e-01 7.92e-01h  1
  18 -1.0125091e+00 2.83e-05 4.85e-04  -6.1 1.06e-01    -  1.00e+00 1.00e+00h  1
  19 -1.0125585e+00 1.21e-05 8.67e-05  -6.8 1.01e-01    -  1.00e+00 1.00e+00h  1
iter    objective    inf_pr   inf_du lg(mu)  ||d||  lg(rg) alpha_du alpha_pr  ls
  20 -1.0125675e+00 4.40e-06 3.01e-05  -7.3 9.80e-02    -  1.00e+00 1.00e+00h  1
  21 -1.0125707e+00 1.75e-06 3.70e-02  -8.0 1.14e-01    -  1.00e+00 9.91e-01h  1
  22 -1.0125714e+00 9.50e-07 4.98e-03  -8.6 1.37e-01    -  1.00e+00 9.94e-01h  1

Number of Iterations....: 22

                                   (scaled)                 (unscaled)
Objective...............:  -1.0125714428815771e+00   -1.0125714428815771e+00
Dual infeasibility......:   4.9826924050513349e-03    4.9826924050513349e-03
Constraint violation....:   3.5941849951814930e-07    9.4995894309168882e-07
Variable bound violation:   9.9621924876114321e-09    9.9621924876114321e-09
Complementarity.........:   8.1796266163751730e-09    8.1796266163751730e-09
Overall NLP error.......:   5.1484604230090588e-07    4.9826924050513349e-03


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

EXIT: Optimal Solution Found.

and plot the solution

plt = plot(direct_sol, solution_label="(direct)", size=(800, 800))
Example block output

Structure of the solution

We first determine visually the structure of the optimal solution which is composed of a bang arc with maximal control, followed by a singular arc, then by a boundary arc and the final arc is with zero control. Note that the switching function vanishes along the singular and boundary arcs.

Interactions with an optimal control solution

Please check state, costate, control and variable to get data from the solution. The functions state, costate and control return functions of time and variable returns a vector. The function time_grid returns the discretized time grid returned by the solver.

t = time_grid(direct_sol)
x = state(direct_sol)
u = control(direct_sol)
p = costate(direct_sol)

H1 = Lift(F1)           # H1(x, p) = p' * F1(x)
φ(t) = H1(x(t), p(t))   # switching function
g(x) = vmax - x[2]      # state constraint v ≤ vmax

u_plot  = plot(t, u,     label = "u(t)")
H1_plot = plot(t, φ,     label = "H₁(x(t), p(t))")
g_plot  = plot(t, g ∘ x, label = "g(x(t))")

plot(u_plot, H1_plot, g_plot, layout=(3,1), size=(500, 500))
Example block output

We are now in position to solve the problem by an indirect shooting method. We first define the four control laws in feedback form and their associated flows. For this we need to compute some Lie derivatives, namely Poisson brackets of Hamiltonians (themselves obtained as lifts to the cotangent bundle of vector fields), or derivatives of functions along a vector field. For instance, the control along the minimal order singular arcs is obtained as the quotient

\[u_s = -\frac{H_{001}}{H_{101}}\]

of length three Poisson brackets:

\[H_{001} = \{H_0,\{H_0,H_1\}\}, \quad H_{101} = \{H_1,\{H_0,H_1\}\}\]

where, for two Hamiltonians $H$ and $G$,

\[\{H,G\} := (\nabla_p H|\nabla_x G) - (\nabla_x H|\nabla_p G).\]

While the Lie derivative of a function $f$ wrt. a vector field $X$ is simply obtained as

\[(X \cdot f)(x) := f'(x) \cdot X(x),\]

and is used to the compute the control along the boundary arc,

\[u_b(x) = -(F_0 \cdot g)(x) / (F_1 \cdot g)(x),\]

as well as the associated multiplier for the order one state constraint on the velocity:

\[\mu(x, p) = H_{01}(x, p) / (F_1 \cdot g)(x).\]

Poisson bracket and Lie derivative

The Poisson bracket $\{H,G\}$ is also given by the Lie derivative of $G$ along the Hamiltonian vector field $X_H = (\nabla_p H, -\nabla_x H)$ of $H$, that is

\[ \{H,G\} = X_H \cdot G\]

which is the reason why we use the @Lie macro to compute Poisson brackets below.

With the help of the differential geometry primitives from CTBase.jl, these expressions are straightforwardly translated into Julia code:

# Controls
u0 = 0                                  # off control
u1 = 1                                  # bang control

H0 = Lift(F0)                           # H0(x, p) = p' * F0(x)
H01  = @Lie { H0, H1 }
H001 = @Lie { H0, H01 }
H101 = @Lie { H1, H01 }
us(x, p) = -H001(x, p) / H101(x, p)     # singular control

ub(x) = -(F0⋅g)(x) / (F1⋅g)(x)          # boundary control
μ(x, p) = H01(x, p) / (F1⋅g)(x)         # multiplier associated to the state constraint g

# Flows
f0 = Flow(ocp, (x, p, tf) -> u0)
f1 = Flow(ocp, (x, p, tf) -> u1)
fs = Flow(ocp, (x, p, tf) -> us(x, p))
fb = Flow(ocp, (x, p, tf) -> ub(x), (x, u, tf) -> g(x), (x, p, tf) -> μ(x, p))

Shooting function

Then, we define the shooting function according to the optimal structure we have determined, that is a concatenation of four arcs.

x0 = [ r0, v0, m0 ] # initial state

function shoot!(s, p0, t1, t2, t3, tf)

    x1, p1 = f1(t0, x0, p0, t1)
    x2, p2 = fs(t1, x1, p1, t2)
    x3, p3 = fb(t2, x2, p2, t3)
    xf, pf = f0(t3, x3, p3, tf)

    s[1] = xf[3] - mf                             # final mass constraint
    s[2:3] = pf[1:2] - [ 1, 0 ]                   # transversality conditions
    s[4] = H1(x1, p1)                             # H1 = H01 = 0
    s[5] = H01(x1, p1)                            # at the entrance of the singular arc
    s[6] = g(x2)                                  # g = 0 when entering the boundary arc
    s[7] = H0(xf, pf)                             # since tf is free

end

Initial guess

To solve the problem by an indirect shooting method, we then need a good initial guess, that is a good approximation of the initial costate, the three switching times and the final time.

η = 1e-3
t13 = t[ abs.(φ.(t)) .≤ η ]
t23 = t[ 0 .≤ (g ∘ x).(t) .≤ η ]
p0 = p(t0)
t1 = min(t13...)
t2 = min(t23...)
t3 = max(t23...)
tf = t[end]

println("p0 = ", p0)
println("t1 = ", t1)
println("t2 = ", t2)
println("t3 = ", t3)
println("tf = ", tf)

# Norm of the shooting function at solution
using LinearAlgebra: norm
s = similar(p0, 7)
shoot!(s, p0, t1, t2, t3, tf)
println("\nNorm of the shooting function: ‖s‖ = ", norm(s), "\n")
p0 = [3.942902223769028, 0.1463010011126008, 0.05413331753584351]
t1 = 0.01817624054194442
t2 = 0.06058746847314807
t3 = 0.09895953183947517
tf = 0.2019582282438269

Norm of the shooting function: ‖s‖ = 3.247075939938923

Indirect shooting

We aggregate the data to define the initial guess vector.

ξ = [ p0 ; t1 ; t2 ; t3 ; tf ] # initial guess
7-element Vector{Float64}:
 3.942902223769028
 0.1463010011126008
 0.05413331753584351
 0.01817624054194442
 0.06058746847314807
 0.09895953183947517
 0.2019582282438269

NonlinearSolve.jl

We first use the NonlinearSolve.jl package to solve the shooting equation. Let us define the problem.

# auxiliary function with aggregated inputs
nle! = (s, ξ, λ) -> shoot!(s, ξ[1:3], ξ[4], ξ[5], ξ[6], ξ[7])

# NLE problem with initial guess
prob = NonlinearProblem(nle!, ξ)

Let us do some benchmarking.

using BenchmarkTools
@benchmark solve(prob; abstol=1e-8, reltol=1e-8, show_trace=Val(false))
BenchmarkTools.Trial: 4 samples with 1 evaluation per sample.
 Range (minmax):  1.482 s  1.592 s   GC (min … max): 23.43% … 13.67%
 Time  (median):     1.534 s               GC (median):    22.02%
 Time  (mean ± σ):   1.536 s ± 56.507 ms   GC (mean ± σ):  20.23% ±  4.61%

  █                                               █      █  
  █▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁█▁▁▁▁▁▁█ ▁
  1.48 s         Histogram: frequency by time        1.59 s <

 Memory estimate: 3.92 GiB, allocs estimate: 1278725.

For small nonlinear systems, it could be faster to use the SimpleNewtonRaphson() descent algorithm.

@benchmark solve(prob, SimpleNewtonRaphson(); abstol=1e-8, reltol=1e-8, show_trace=Val(false))
BenchmarkTools.Trial: 4 samples with 1 evaluation per sample.
 Range (minmax):  1.369 s  1.426 s   GC (min … max): 23.38% … 11.19%
 Time  (median):     1.396 s               GC (median):    23.09%
 Time  (mean ± σ):   1.397 s ± 28.604 ms   GC (mean ± σ):  20.33% ±  6.14%

  █     █                                        █        █  
  █▁▁▁▁▁█▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁█▁▁▁▁▁▁▁▁█ ▁
  1.37 s         Histogram: frequency by time        1.43 s <

 Memory estimate: 3.62 GiB, allocs estimate: 1081343.

Now, let us solve the problem and retrieve the initial costate solution.

# resolution of S(ξ) = 0
indirect_sol = solve(prob; abstol=1e-8, reltol=1e-8, show_trace=Val(true))

# we retrieve the costate solution together with the times
p0 = indirect_sol.u[1:3]
t1 = indirect_sol.u[4]
t2 = indirect_sol.u[5]
t3 = indirect_sol.u[6]
tf = indirect_sol.u[7]

println("")
println("p0 = ", p0)
println("t1 = ", t1)
println("t2 = ", t2)
println("t3 = ", t3)
println("tf = ", tf)

# Norm of the shooting function at solution
s = similar(p0, 7)
shoot!(s, p0, t1, t2, t3, tf)
println("\nNorm of the shooting function: ‖s‖ = ", norm(s), "\n")

Algorithm: NewtonRaphson(
    descent = NewtonDescent(),
    autodiff = AutoForwardDiff(),
    vjp_autodiff = AutoReverseDiff(
        compile = false
    ),
    jvp_autodiff = AutoForwardDiff(),
    concrete_jac = Val{false}()
)

----    	-------------       	-----------
Iter    	f(u) inf-norm       	Step 2-norm
----    	-------------       	-----------
0       	3.24457169e+00      	NaN
1       	2.37667368e-01      	1.99209184e-01
2       	3.48822712e-02      	2.01514150e-01
3       	2.77063016e-04      	7.88857353e-04
4       	4.22201085e-09      	6.46398261e-05
Final   	4.22201085e-09
----------------------

p0 = [3.9457646578192955, 0.1503955962287936, 0.053712712941242795]
t1 = 0.02350968401827331
t2 = 0.05973738093914581
t3 = 0.10157134843292695
tf = 0.20204744062225508

Norm of the shooting function: ‖s‖ = 4.2244265609382325e-9

MINPACK.jl

Instead of the NonlinearSolve.jl package we can use the MINPACK.jl package to solve the shooting equation. To compute the Jacobian of the shooting function we use the DifferentiationInterface.jl package with ForwardDiff.jl backend.

using DifferentiationInterface
import ForwardDiff
backend = AutoForwardDiff()

Let us define the problem to solve.

# auxiliary function with aggregated inputs
nle!  = ( s, ξ) -> shoot!(s, ξ[1:3], ξ[4], ξ[5], ξ[6], ξ[7])

# Jacobian of the (auxiliary) shooting function
jnle! = (js, ξ) -> jacobian!(nle!, similar(ξ), js, backend, ξ)

We are now in position to solve the problem with the hybrj solver from MINPACK.jl through the fsolve function, providing the Jacobian. Let us do some benchmarking.

@benchmark fsolve(nle!, jnle!, ξ; show_trace=false) # initial guess given to the solver
BenchmarkTools.Trial: 9 samples with 1 evaluation per sample.
 Range (minmax):  522.270 ms697.753 ms   GC (min … max): 15.31% … 7.44%
 Time  (median):     601.884 ms                GC (median):    22.47%
 Time  (mean ± σ):   605.363 ms ±  46.553 ms   GC (mean ± σ):  20.66% ± 5.74%

  ▁                  ▁ ▁     █    ▁ ▁▁                        ▁  
  █▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁█▁█▁▁▁▁▁█▁▁█▁██▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁█ ▁
  522 ms           Histogram: frequency by time          698 ms <

 Memory estimate: 1.50 GiB, allocs estimate: 1184652.

We can also use the preparation step of DifferentiationInterface.jl.

extras = prepare_jacobian(nle!, similar(ξ), backend, ξ)
jnle_prepared!(js, ξ) = jacobian!(nle!, similar(ξ), js, backend, ξ, extras)
@benchmark fsolve(nle!, jnle_prepared!, ξ; show_trace=false)
BenchmarkTools.Trial: 9 samples with 1 evaluation per sample.
 Range (minmax):  525.366 ms609.758 ms   GC (min … max): 23.59% … 6.80%
 Time  (median):     554.930 ms                GC (median):    25.17%
 Time  (mean ± σ):   559.154 ms ±  24.922 ms   GC (mean ± σ):  22.06% ± 6.46%

  █    █           ██  █      █   █                       █  
  █▁▁▁▁█▁▁▁▁▁▁▁▁▁▁▁██▁█▁▁█▁▁▁▁▁▁█▁▁▁█▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁▁█ ▁
  525 ms           Histogram: frequency by time          610 ms <

 Memory estimate: 1.28 GiB, allocs estimate: 1815902.

Now, let us solve the problem and retrieve the initial costate solution.

# resolution of S(ξ) = 0
indirect_sol = fsolve(nle!, jnle!, ξ, show_trace=true)

# we retrieve the costate solution together with the times
p0 = indirect_sol.x[1:3]
t1 = indirect_sol.x[4]
t2 = indirect_sol.x[5]
t3 = indirect_sol.x[6]
tf = indirect_sol.x[7]

println("")
println("p0 = ", p0)
println("t1 = ", t1)
println("t2 = ", t2)
println("t3 = ", t3)
println("tf = ", tf)

# Norm of the shooting function at solution
s = similar(p0, 7)
shoot!(s, p0, t1, t2, t3, tf)
println("\nNorm of the shooting function: ‖s‖ = ", norm(s), "\n")
Iter     f(x) inf-norm    Step 2-norm      Step time
------   --------------   --------------   --------------
     1     3.244572e+00     0.000000e+00         0.252568
     2     2.376674e-01     3.968430e-02         0.446647
     3     1.323103e-01     1.239656e-02         0.017537
     4     6.581010e-02     9.419374e-03         0.026899
     5     2.924902e-03     1.078664e-04         0.025393
     6     2.649943e-03     6.305981e-06         0.026598
     7     2.181309e-04     1.692478e-06         0.026368
     8     8.344768e-05     5.907443e-08         0.026390
     9     2.221351e-05     9.554355e-11         0.018136
    10     3.390196e-06     2.715857e-11         0.040055
    11     9.798531e-08     4.035506e-13         0.025325
    12     1.005096e-08     1.846535e-16         0.025315

p0 = [3.945764657131959, 0.15039559627320845, 0.05371271297238469]
t1 = 0.02350968397201485
t2 = 0.05973738106803653
t3 = 0.10157134845062943
tf = 0.20204744056815385

Norm of the shooting function: ‖s‖ = 1.0077450807057513e-8

Plot of the solution

We plot the solution of the indirect solution (in red) over the solution of the direct method (in blue).

f = f1 * (t1, fs) * (t2, fb) * (t3, f0) # concatenation of the flows
flow_sol = f((t0, tf), x0, p0)          # compute the solution: state, costate, control...

plot!(plt, flow_sol, solution_label="(indirect)")
Example block output

References

  • 1R.H. Goddard. A Method of Reaching Extreme Altitudes, volume 71(2) of Smithsonian Miscellaneous Collections. Smithsonian institution, City of Washington, 1919.
  • 2H. Seywald and E.M. Cliff. Goddard problem in presence of a dynamic pressure limit. Journal of Guidance, Control, and Dynamics, 16(4):776–781, 1993.