Nonlinear Optimization Solution
Nonlinear Optimization solution
uses generateNLJuMPcodes
to solve the Problems. The arguments are the same as generateJuMPcodes
OptControl.generateNLJuMPcodes
— FunctiongenerateNLJuMPcodes(L, F, state, u, tspan, t0, tf) -> Any
generateNLJuMPcodes(L, F, state, u, tspan, t0, tf, Φ) -> Any
generateNLJuMPcodes(
L,
F,
state,
u,
tspan,
t0,
tf,
Φ,
tf_constraint;
state_ub,
state_lb,
u_ub,
u_lb,
N,
discretization,
model_name,
writeFilePath,
tolerance
) -> Any
generateNLJuMPcodes
is used to solve Nonlinear control problem
A general form of optimal control problem:
\[min (\Phi(\boldsymbol{x}(t_f),t_f)+\int_{t_0}^{t_f} L[\boldsymbol{x}(t),\boldsymbol{u}(t),t]dt) \\ s.t. \hspace{0.4cm} \dot{\boldsymbol{x}} =f[\boldsymbol{x}(t),\boldsymbol{u}(t),t]\]
args:
L
: above L
\[L[\boldsymbol{x}(t),\boldsymbol{u}(t),t]\]
F
: above f
\[f[\boldsymbol{x}(t),\boldsymbol{u}(t),t]\]
state
: states(variable) inF
u
: control variable u inF
tspan
: time fieldt0
: initial value attspan[1]
, length oft0
must be equal to length of statetf
: final(end) value attspan[2]
, length oftf
must be equal to length of stateΦ
: above Φ, default: nothing
\[\Phi(\boldsymbol{x}(t_f),t_f)\]
tf_constraint
: special requirements for end time, default: nothing. Example:
\[x_{1}+x_{2}=0\]
state_ub
: state's upper limit, length ofstate_ub
must be equal to length of state, default: nothingstate_lb
: state's lower limit, length ofstate_lb
must be equal to length of state, default: nothingu_ub
: u's upper limit, length ofu_ub
must be equal to length of u, default: nothingu_lb
: u's lower limit, length ofu_lb
must be equal to length of u, default: nothingN
: Number of discrete, default: 1000
\[dt = (endTime - startTime) / N\]
discretization
: discretization methods, default: "trapezoidal"model_name
: name of JuMP model, default: "model"writeFilePath
: path of generated code, default: nothingtolerance
: acceptable_tol of Ipopt, default: 1.0e-6
Example:
using OptControl,ModelingToolkit
@variables t u x
f = exp(x) + u
L = 0.5 * u^2
t0 = [1.0]
tf = [0.0]
tspan = (0.0, 2.0)
N = 100
sol = generateNLJuMPcodes(L, f, x, u, tspan, t0, tf; N=N)
Let's see some examples. If you are not familar with ModelingToolkit.jl, you can see the Symbolics.jl document (which ModelingToolkit.jl based on) and try some tests about symbolic computation
using OptControl,ModelingToolkit
@variables x[1:2] y[1:2]
print(scalarize(rand(1:10,2,2)*x+rand(1:10,2,2)*y))
Num[3x[2] + 4y[1] + 6y[2] + x[1], 2x[1] + 4x[2] + 6y[1] + 6y[2]]
PS: scalarize is from Symbolics.jl
If you need, pass a name to writeFilePath
and do some changes in script.
Example 1
To solve:
\[min \int_{0}^{2} u^2dt \newline s.t. ~~~~~ \dot{x} =e^x+ u \newline x(0) = 1, x(2)=0\]
Just define variables and build functions. Call generateNLJuMPcodes
and get the results.
using OptControl, ModelingToolkit, Test
@variables t u x
f = exp(x) + u
L = u^2
t0 = [1.0]
tf = [0.0]
tspan = (0.0, 2.0)
N = 100
sol = generateNLJuMPcodes(L, f, x, u, tspan, t0, tf; N=N)
@test isapprox.(0.0, sol[1][end, :], atol=1.0e-5) == [true for i in 1:length(sol[1][end, :])]
Test Passed
Example 2
To solve:
\[min \int_{0}^{2} u^2dt \newline s.t. ~~~~~ \dot{\boldsymbol{x}} =\begin{bmatrix}exp&cos \newline sin&1\end{bmatrix}\boldsymbol{x}+ \begin{bmatrix}0 \newline 1 \end{bmatrix}u \newline \boldsymbol{x}(0) = \begin{bmatrix} 1 \newline 1 \end{bmatrix}, \boldsymbol{x}(1)=\begin{bmatrix} 0 \newline 0 \end{bmatrix}\]
using OptControl, ModelingToolkit, Test
@variables t u x[1:2]
f = [exp(x[1]) + cos(x[2]), sin(x[1]) + x[2]] + [1, 0] * u
L = u^2
t0 = [1.0, 1.0]
tf = [0.0, 0.0]
tspan = (0.0, 2.0)
N = 100
sol = generateNLJuMPcodes(L, f, x, u, tspan, t0, tf; N=N)
@test isapprox.(0.0, sol[1][end, :], atol=1.0e-5) == [true for i in 1:length(sol[1][end, :])]
Test Passed