Bioreactor

This problem models a coupled photobioreactor–digester system for methane production. The system consists of three state variables: the algae concentration $y(t)$, the substrate concentration $s(t)$, and the biomass concentration $b(t)$. The control variable $u(t)$ represents the input flow rate between the two units. The dynamics include algal growth driven by light, substrate consumption, and biomass evolution. The aim is to maximise methane production over a fixed time horizon under biological and operational constraints.

Problem formulation

We minimise

\[\min_{y,\,s,\,b,\,u} J(y,s,b,u) = - \frac{1}{\beta+c} \int_0^T \mu_2(s(t))\, b(t)\, dt,\]

subject to the dynamics

\[\dot{y}(t) = \frac{\mu(t)\, y(t)}{1+y(t)} - (r+u(t))\, y(t),\]

\[\dot{s}(t) = -\mu_2(s(t))\, b(t) + u(t)\, \beta\big(\gamma y(t) - s(t)\big),\]

\[\dot{b}(t) = \big(\mu_2(s(t)) - u(t)\, \beta\big)\, b(t),\]

with

  • Light model: $\mu(t) = \mu_{\text{bar}} \,\max(0,\sin(\tau(t)))^2$, where $\tau(t)$ encodes a periodic day–night cycle,
  • Growth function (Monod law): $\mu_2(s) = \mu_2^m \,\dfrac{s}{K_s + s}$.

Constraints

  • Control bounds: $0 \leq u(t) \leq 1$,
  • State bounds: $y(t) \geq 0,\; s(t) \geq 0,\; b(t) \geq 10^{-3}$,
  • Initial conditions: $0.05 \leq y(0) \leq 0.25$, $0.5 \leq s(0) \leq 5$, $0.5 \leq b(0) \leq 3$.

The horizon is fixed to $T = 200$ (rescaled units), corresponding to several day–night periods.

Qualitative behaviour

The optimal solution exploits the periodic structure of the problem: during illuminated phases, algal growth is favoured, while in dark phases methane production becomes predominant. The control exhibits a near bang–bang structure, alternating between minimal and maximal input flows, with possible singular arcs when substrate and biomass reach balanced levels. The constraints on positivity of the states are typically active for $b(t)$ at the beginning of the process.

References

  • Bayen, T., Mairet, F., Martinon, P., & Sebbah, M. (2014). Analysis of a periodic optimal control problem connected to microalgae anaerobic digestion. Optimal Control Applications and Methods. DOI:10.1002/oca.2127
  • Bayen, T., Mairet, F., Martinon, P., & Sebbah, M. (2013). Optimising the anaerobic digestion of microalgae in a coupled process. 13th European Control Conference.
  • Barbosa, M.J., & Wijffels, R.H. (2010). An Outlook on Microalgal Biofuels. Science, 329, 796–799.
  • Betts, J.T. (2001). Practical methods for optimal control using nonlinear programming. SIAM.
  • BOCOP repository: https://github.com/control-toolbox/bocop/tree/main/bocop

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(:bioreactor)
Example block output

Solve the problem

OptimalControl model

Import the OptimalControl model and solve it.

# import DOCP model
docp = bioreactor(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............................:     2404
                     variables with only lower bounds:     1803
                variables with lower and upper bounds:      601
                     variables with only upper bounds:        0
Total number of equality constraints.................:     1800
Total number of inequality constraints...............:        3
        inequality constraints with only lower bounds:        0
   inequality constraints with lower and upper bounds:        3
        inequality constraints with only upper bounds:        0


Number of Iterations....: 169

                                   (scaled)                 (unscaled)
Objective...............:  -1.2925452053485811e+01   -1.2925452053485811e+01
Dual infeasibility......:   1.7430226191149280e-09    1.7430226191149280e-09
Constraint violation....:   2.7871038810189930e-11    2.7871038810189930e-11
Variable bound violation:   9.9794036661269432e-09    9.9794036661269432e-09
Complementarity.........:   5.6482249599314883e-11    5.6482249599314883e-11
Overall NLP error.......:   1.7430226191149280e-09    1.7430226191149280e-09


Number of objective function evaluations             = 523
Number of objective gradient evaluations             = 170
Number of equality constraint evaluations            = 523
Number of inequality constraint evaluations          = 523
Number of equality constraint Jacobian evaluations   = 170
Number of inequality constraint Jacobian evaluations = 170
Number of Lagrangian Hessian evaluations             = 169
Total seconds in IPOPT                               = 2.846

EXIT: Optimal Solution Found.

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

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

JuMP model

Import the JuMP model and solve it.

# import model
nlp_jp = bioreactor(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............................:     2404
                     variables with only lower bounds:     1803
                variables with lower and upper bounds:      601
                     variables with only upper bounds:        0
Total number of equality constraints.................:     1800
Total number of inequality constraints...............:        3
        inequality constraints with only lower bounds:        0
   inequality constraints with lower and upper bounds:        3
        inequality constraints with only upper bounds:        0


Number of Iterations....: 199

                                   (scaled)                 (unscaled)
Objective...............:  -1.2925452054139990e+01   -1.2925452054139990e+01
Dual infeasibility......:   2.1331552219061578e-09    2.1331552219061578e-09
Constraint violation....:   3.3919533848347783e-12    3.3919533848347783e-12
Variable bound violation:   9.9843364557546566e-09    9.9843364557546566e-09
Complementarity.........:   2.3792290663209013e-11    2.3792290663209013e-11
Overall NLP error.......:   2.1331552219061578e-09    2.1331552219061578e-09


Number of objective function evaluations             = 739
Number of objective gradient evaluations             = 200
Number of equality constraint evaluations            = 739
Number of inequality constraint evaluations          = 739
Number of equality constraint Jacobian evaluations   = 200
Number of inequality constraint Jacobian evaluations = 200
Number of Lagrangian Hessian evaluations             = 199
Total seconds in IPOPT                               = 5.759

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_order169-12.9255
2JuMPLOCALLY_SOLVED199-12.9255

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(:bioreactor, docp, nlp_oc_sol, nlp_jp)
┌─ bioreactor
│
├─  Number of iterations
│
│     OptimalControl : 169
│     JuMP           : 199
│
├─  State y (L2 norm)
│
│     Absolute error : 5.71160038272571e-7
│     Relative error : 6.992660113504655e-9
│
├─  State s (L2 norm)
│
│     Absolute error : 2.82699948432636e-7
│     Relative error : 1.6118984742550558e-8
│
├─  State b (L2 norm)
│
│     Absolute error : 6.229457079560391e-7
│     Relative error : 7.511859329693764e-9
│
├─  Control u (L2 norm)
│
│     Absolute error : 2.0170976125584572e-6
│     Relative error : 2.4263528710751327e-6
│
├─  objective
│
│     Absolute error : 6.541789332459302e-10
│     Relative error : 5.0611686966333444e-11
│
└─

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[:bioreactor][:state_name])
m = control_dimension(ocp_sol) # or length(metadata[:bioreactor][: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(:bioreactor, nlp_jp)     # t0, ..., tN = tf
x = state(:bioreactor, nlp_jp)         # function of time
u = control(:bioreactor, nlp_jp)       # function of time
p = costate(:bioreactor, 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