Brachistochrone
Description of the problem
The Brachistochrone problem is a classical benchmark in the history of calculus of variations and optimal control. It seeks the shape of a curve (or wire) connecting two points $A$ and $B$ within a vertical plane, such that a bead sliding frictionlessly under the influence of uniform gravity travels from $A$ to $B$ in the shortest possible time. Originating from the challenge posed by Johann Bernoulli in 1696 Bernoulli 1696, it marks the birth of modern optimal control theory. The problem involves nonlinear dynamics where the state includes position and velocity, and the control is the path's angle, with the final time $t_f$ being a decision variable to be minimised.
Mathematical formulation
The problem can be stated as a free final time optimal control problem:
\[\begin{aligned} \min_{u(\cdot), t_f} \quad & J = t_f = \int_0^{t_f} 1 \,\mathrm{d}t \\[1em] \text{s.t.} \quad & \dot{x}(t) = v(t) \sin u(t), \\[0.5em] & \dot{y}(t) = v(t) \cos u(t), \\[0.5em] & \dot{v}(t) = g \cos u(t), \\[0.5em] & x(0) = x_0, \; y(0) = y_0, \; v(0) = v_0, \\[0.5em] & x(t_f) = x_f, \; y(t_f) = y_f, \end{aligned}\]
where $u(t)$ represents the angle of the tangent to the curve with respect to the vertical axis, and $g$ is the gravitational acceleration.
Qualitative behaviour
This problem is a classic example of minimum time control with nonlinear dynamics. Unlike the shortest path problem (a straight line), the optimal trajectory balances the need to minimize path length with the need to maximize velocity.
The analytical solution to this problem is a cycloid—the curve traced by a point on the rim of a circular wheel as the wheel rolls along a straight line without slipping. Specifically:
- Energy Conservation: Since the system is conservative and frictionless, the speed at any height $h$ is given by $v = \sqrt{2gh}$ (assuming start from rest). This implies that maximizing vertical drop early in the trajectory increases velocity for the remainder of the path.
- Concavity: The optimal curve is concave up. It starts steeply (vertical tangent if $v_0=0$) to gain speed quickly, then flattens out near the bottom before potentially rising again to reach the target.
- Singularity: At the initial point, if $v_0=0$, the control is theoretically singular as the bead has no velocity to direct. Numerical solvers often require a small non-zero initial velocity or a robust initial guess to handle this start.
The control $u(t)$ is continuous and varies smoothly to trace the cycloidal arc.
Characteristics
- Nonlinear dynamics involving trigonometric functions of the control.
- Free final time problem (time-optimal).
- Analytical solution is a Cycloid.
- Historically significant as the first problem solved using techniques that evolved into the Calculus of Variations.
References
- Bernoulli, J. (1696). Problema novum ad cujus solutionem Mathematici invitantur.
Acta Eruditorum. The original publication posing the challenge to the mathematicians of Europe, solved by Newton, Leibniz, L'Hôpital, and the Bernoulli brothers.
- Sussmann, H. J., & Willems, J. C. (1997). 300 years of optimal control: from the brachystochrone to the maximum principle.
IEEE Control Systems Magazine. doi.org/10.1109/37.588108 A comprehensive historical review linking the classical Brachistochrone problem to modern optimal control theory and Pontryagin's Maximum Principle.
- Dymos Examples: Brachistochrone. OpenMDAO/Dymos
The numerical implementation serving as the basis for the current benchmark formulation.
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(:brachistochrone)[:grid_size]500The default values of the parameters are:
metadata(:brachistochrone)[:parameters]Parameter = Value
------------------
g = 9.8066e+00
t0 = 0.0000e+00
x0 = 0.0000e+00
y0 = 1.0000e+01
v0 = 0.0000e+00
xf = 1.0000e+01
yf = 5.0000e+00
u_min = -3.1416e+00
u_max = 3.1416e+00Initial 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.
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
endplot_initial_guess(:brachistochrone)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=:brachistochrone,
Grid_Size=metadata(:brachistochrone)[:grid_size],
Variables=get_nvar(nlp_model(brachistochrone(OptimalControlBackend()))),
Constraints=get_ncon(nlp_model(brachistochrone(OptimalControlBackend()))),
))| Row | Problem | Grid_Size | Variables | Constraints |
|---|---|---|---|---|
| Symbol | Int64 | Int64 | Int64 | |
| 1 | brachistochrone | 500 | 2005 | 1505 |
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 = brachistochrone(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: 0
variables with lower and upper bounds: 502
variables with only upper bounds: 0
Total number of equality constraints.................: 1505
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....: 18
(scaled) (unscaled)
Objective...............: 1.8029320066529224e+00 1.8029320066529224e+00
Dual infeasibility......: 7.9491352468077393e-10 7.9491352468077393e-10
Constraint violation....: 6.4120779796894567e-09 6.4120779796894567e-09
Variable bound violation: 2.6279822762376170e-09 2.6279822762376170e-09
Complementarity.........: 1.6184387283906659e-09 1.6184387283906659e-09
Overall NLP error.......: 6.4120779796894567e-09 6.4120779796894567e-09
Number of objective function evaluations = 20
Number of objective gradient evaluations = 19
Number of equality constraint evaluations = 20
Number of inequality constraint evaluations = 0
Number of equality constraint Jacobian evaluations = 19
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations = 18
Total seconds in IPOPT = 1.462
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 = brachistochrone(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: 0
variables with lower and upper bounds: 502
variables with only upper bounds: 0
Total number of equality constraints.................: 1505
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....: 18
(scaled) (unscaled)
Objective...............: 1.8029320066529222e+00 1.8029320066529222e+00
Dual infeasibility......: 7.9491496797070594e-10 7.9491496797070594e-10
Constraint violation....: 6.4120770915110370e-09 6.4120770915110370e-09
Variable bound violation: 2.6279822762376170e-09 2.6279822762376170e-09
Complementarity.........: 1.6184387282039994e-09 1.6184387282039994e-09
Overall NLP error.......: 6.4120770915110370e-09 6.4120770915110370e-09
Number of objective function evaluations = 20
Number of objective gradient evaluations = 19
Number of equality constraint evaluations = 20
Number of inequality constraint evaluations = 0
Number of equality constraint Jacobian evaluations = 19
Number of inequality constraint Jacobian evaluations = 0
Number of Lagrangian Hessian evaluations = 18
Total seconds in IPOPT = 0.143
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),
))| Row | Model | Flag | Iterations | Objective |
|---|---|---|---|---|
| Symbol | Any | Int64 | Float64 | |
| 1 | OptimalControl | first_order | 18 | 1.80293 |
| 2 | JuMP | LOCALLY_SOLVED | 18 | 1.80293 |
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.
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
endprint_numerical_comparisons(:brachistochrone, docp, nlp_oc_sol, nlp_jp)┌─ brachistochrone
│
├─ Number of Iterations
│ OptimalControl : 18 JuMP : 18
├─ States (L2 Norms)
│ px Abs: 3.891e-15 Rel: 6.568e-16
│ py Abs: 2.083e-15 Rel: 2.088e-16
│ v Abs: 9.599e-15 Rel: 9.647e-16
├─ Controls (L2 Norms)
│ u Abs: 2.045e-15 Rel: 1.508e-15
├─ Variables
│ tf Abs: 2.220e-16 Rel: 1.232e-16
├─ Objective
│ Abs: 2.220e-16 Rel: 1.232e-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