Skip to content

Commit ae36b96

Browse files
init struct and constructors (constant init, warm start) (#48)
1 parent 815d5a6 commit ae36b96

File tree

8 files changed

+145
-46
lines changed

8 files changed

+145
-46
lines changed

src/CTDirect.jl

Lines changed: 2 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -13,6 +13,7 @@ const __display = CTBase.__display
1313
const matrix2vec = CTBase.matrix2vec
1414

1515
# includes
16+
include("init.jl")
1617
include("utils.jl")
1718
include("problem.jl")
1819
include("solution.jl")
@@ -21,5 +22,6 @@ include("solve.jl")
2122
# export functions only for user
2223
export solve
2324
export is_solvable
25+
export OptimalControlInit
2426

2527
end

src/init.jl

Lines changed: 52 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,52 @@
1+
# NB. to be moved up to CTBase later
2+
3+
# init struct (to be passed to solve)
4+
# functions of time for state and control variables
5+
# values for optimization variables (NB. free t0/tf in particular !)
6+
# later for indirect methods add costate and structure information
7+
8+
# constructors
9+
# - default: keep empty for method-specific handling
10+
# - from constant vector (dim x + dim u + dim v)
11+
# - from existing solution
12+
mutable struct OptimalControlInit
13+
14+
state_dimension
15+
control_dimension
16+
variable_dimension
17+
state_init
18+
control_init
19+
variable_init
20+
info
21+
22+
function OptimalControlInit()
23+
init = new()
24+
init.info = :undefined
25+
return init
26+
end
27+
28+
function OptimalControlInit(x_init, u_init, v_init=Real[])
29+
init = new()
30+
init.info = :constant
31+
init.state_dimension = length(x_init)
32+
init.control_dimension = length(u_init)
33+
init.variable_dimension = length(v_init)
34+
#println("Init dims x u v: ",init.state_dimension," ",init.control_dimension," ",init.variable_dimension)
35+
init.state_init = t -> x_init
36+
init.control_init = t -> u_init
37+
init.variable_init = v_init
38+
return init
39+
end
40+
41+
function OptimalControlInit(sol::OptimalControlSolution)
42+
init = new()
43+
init.info = :solution
44+
init.state_dimension = sol.state_dimension
45+
init.control_dimension = sol.control_dimension
46+
init.variable_dimension = sol.variable_dimension
47+
init.state_init = t-> sol.state(t)[1:sol.state_dimension]
48+
init.control_init = sol.control
49+
init.variable_init = sol.infos[:variable] #+++ need proper variable field in ocp solution !
50+
return init
51+
end
52+
end

src/problem.jl

Lines changed: 9 additions & 4 deletions
Original file line numberDiff line numberDiff line change
@@ -54,10 +54,12 @@ mutable struct CTDirect_data
5454

5555
## NLP
5656
# NLP problem
57-
dim_NLP_state # 'augmented state'
57+
dim_NLP_state # 'augmented state' with possible additional component for lagrange cost
5858
dim_NLP_constraints
5959
dim_NLP_variables
6060
dim_NLP_steps
61+
62+
# initialization
6163
NLP_init
6264

6365
# NLP solution
@@ -68,7 +70,8 @@ mutable struct CTDirect_data
6870
NLP_iterations
6971
NLP_stats # remove later ? type is https://juliasmoothoptimizers.github.io/SolverCore.jl/stable/reference/#SolverCore.GenericExecutionStats
7072

71-
function CTDirect_data(ocp::OptimalControlModel, N::Integer, init=nothing)
73+
74+
function CTDirect_data(ocp::OptimalControlModel, N::Integer, init::OptimalControlInit)
7275

7376
ctd = new()
7477

@@ -224,12 +227,14 @@ end
224227
function variables_bounds(ctd)
225228

226229
N = ctd.dim_NLP_steps
227-
l_var = -Inf*ones(ctd.dim_NLP_variables)
228-
u_var = Inf*ones(ctd.dim_NLP_variables)
230+
l_var = -Inf * ones(ctd.dim_NLP_variables)
231+
u_var = Inf * ones(ctd.dim_NLP_variables)
229232

230233
# NLP variables layout: [X0, X1 .. XN, U0, U1 .. UN, V]
231234
# NB. keep offset for each block since blocks are optional !
232235

236+
# +++ we could use the setters here (build local vectors then call setter ?!)
237+
233238
# state box
234239
offset = 0
235240
if ctd.has_state_box

src/solution.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -30,7 +30,7 @@ function _OptimalControlSolution(ocp, ipopt_solution, ctd)
3030
sol.state = t -> x(t)
3131
sol.costate = t -> p(t)
3232
sol.control = t -> u(t)
33-
#sol.variable = v +++no field variable?
33+
#sol.variable = v +++need proper field variable
3434
sol.infos[:variable] = v
3535
sol.objective = ctd.NLP_objective
3636
sol.iterations = ctd.NLP_iterations

src/solve.jl

Lines changed: 1 addition & 1 deletion
Original file line numberDiff line numberDiff line change
@@ -32,7 +32,7 @@ function solve(ocp::OptimalControlModel,
3232
print_level::Integer=__print_level_ipopt(),
3333
mu_strategy::String=__mu_strategy_ipopt(),
3434
display::Bool=__display(),
35-
init=nothing,
35+
init::OptimalControlInit=OptimalControlInit(),
3636
kwargs...)
3737

3838
# get full description from partial

src/utils.jl

Lines changed: 41 additions & 32 deletions
Original file line numberDiff line numberDiff line change
@@ -77,72 +77,81 @@ end
7777
function get_final_time(xu, ctd)
7878
if ctd.has_free_final_time
7979
v = get_variable(xu, ctd)
80-
#println("v: ",v)
81-
#println("ctd.final_time: ",ctd.final_time)
82-
#println("v[ctd.final_time]: ", v[ctd.final_time])
8380
return v[ctd.final_time]
8481
else
8582
return ctd.final_time
8683
end
8784
end
8885

86+
#=
87+
function get_time_step(xu, ctd, i)
88+
N = ctd.dim_NLP_steps
89+
@assert i <= N "trying to get time step for i > N"
90+
t0 = get_initial_time(xu, ctd)
91+
tf = get_final_time(xu, ctd)
92+
return t0 + i * (tf - t0) / N
93+
end
94+
=#
8995

9096
## Initialization for the NLP problem
9197

9298
function set_state_at_time_step!(xu, x_init, ctd, i)
9399
nx = ctd.dim_NLP_state
100+
n = ctd.state_dimension
94101
N = ctd.dim_NLP_steps
95102
@assert i <= N "trying to set init for x(t_i) with i > N"
96-
xu[1+i*nx:(i+1)*nx] = x_init[1:nx]
103+
# NB. only set first n components of state variable (nx = n+1 for lagrange cost)
104+
if n == 1
105+
xu[i*nx + 1] = x_init[]
106+
else
107+
xu[i*nx + 1 : i*nx + n] = x_init
108+
end
97109
end
98110

99111
function set_control_at_time_step!(xu, u_init, ctd, i)
100112
nx = ctd.dim_NLP_state
101113
m = ctd.control_dimension
102114
N = ctd.dim_NLP_steps
103115
@assert i <= N "trying to set init for u(t_i) with i > N"
104-
xu[1+(N+1)*nx+i*m:m+(N+1)*nx+i*m] = u_init[1:m]
116+
offset = (N+1)*nx
117+
if m == 1
118+
xu[offset + i*m + 1] = u_init[]
119+
else
120+
xu[offset + i*m + 1 : offset + i*m + m] = u_init
121+
end
105122
end
106123

107124
function set_variable!(xu, v_init, ctd)
108-
xu[end-ctd.variable_dimension+1:end] = v_init[1:ctd.variable_dimension]
125+
if ctd.variable_dimension == 1
126+
xu[end] = v_init[]
127+
else
128+
xu[end-ctd.variable_dimension+1 : end] = v_init
129+
end
109130
end
110131

111132
function initial_guess(ctd)
112133

113-
N = ctd.dim_NLP_steps
114-
init = ctd.NLP_init
134+
# default initialization
135+
# note: internal variables (lagrange cost, k_i for RK schemes) will keep these default values
136+
xu0 = 0.1 * ones(ctd.dim_NLP_variables)
115137

116-
if init === nothing
117-
# default initialization
118-
xu0 = 0.1*ones(ctd.dim_NLP_variables)
119-
else
120-
# split state / control init values
121-
if length(init) != (ctd.state_dimension + ctd.control_dimension + ctd.variable_dimension)
122-
error("vector for initialization should be of size dim_x + dim_u + dim_v ie:",ctd.state_dimension+ctd.control_dimension+ctd.variable_dimension)
123-
end
124-
x_init = zeros(ctd.dim_NLP_state)
125-
x_init[1:ctd.state_dimension] = init[1:ctd.state_dimension]
126-
u_init = zeros(ctd.control_dimension)
127-
u_init[1:ctd.control_dimension] = init[ctd.state_dimension+1:ctd.state_dimension+ctd.control_dimension]
128-
129-
# mayer -> lagrange additional state
130-
if ctd.has_lagrange_cost
131-
x_init[ctd.dim_NLP_state] = 0.1
132-
end
138+
init = ctd.NLP_init
139+
if init.info != :undefined
140+
N = ctd.dim_NLP_steps
141+
t0 = get_initial_time(xu0, ctd)
142+
tf = get_final_time(xu0, ctd)
143+
h = (tf - t0) / N
133144

134-
# set constant initialization for state / control variables
135-
xu0 = zeros(ctd.dim_NLP_variables)
145+
# set state / control variables
136146
for i in 0:N
137-
set_state_at_time_step!(xu0, x_init, ctd, i)
138-
set_control_at_time_step!(xu0, u_init, ctd, i)
147+
ti = t0 + i * h
148+
set_state_at_time_step!(xu0, init.state_init(ti), ctd, i)
149+
set_control_at_time_step!(xu0, init.control_init(ti), ctd, i)
139150
end
140151

141152
# set variables
142153
if (ctd.variable_dimension > 0)
143-
v_init = zeros(ctd.variable_dimension)
144-
v_init[1:ctd.variable_dimension] = init[ctd.state_dimension+ctd.control_dimension+1:ctd.state_dimension+ctd.control_dimension+ctd.variable_dimension]
145-
set_variable!(xu0, v_init, ctd)
154+
set_variable!(xu0, init.variable_init, ctd)
146155
end
147156
end
148157

test/suite/goddard_all_constraints.jl

Lines changed: 22 additions & 6 deletions
Original file line numberDiff line numberDiff line change
@@ -50,10 +50,26 @@ function F1(x)
5050
return [ 0, Tmax/m, -b*Tmax ]
5151
end
5252
dynamics!(ocp, (x, u, v) -> F0(x) + u*F1(x) )
53-
init_constant = [1.05, 0.1, 0.8, 0.5, 0.1]
54-
sol = solve(ocp, grid_size=30, print_level=0, tol=1e-8, init=init_constant)
55-
#println(sol.objective)
56-
#plot(sol)
53+
sol1 = solve(ocp, grid_size=30, print_level=0, tol=1e-8)
5754
@testset verbose = true showtiming = true ":goddard :max_rf :all_constraints_types" begin
58-
@test sol.objective 1.0125 rtol=1e-2
59-
end
55+
@test sol1.objective 1.0125 rtol=1e-2
56+
end
57+
58+
59+
# with constant initial guess
60+
x_init = [1.05, 0.1, 0.8]
61+
u_init = 0.5
62+
v_init = 0.1
63+
init_constant = OptimalControlInit(x_init, u_init, v_init)
64+
sol2 = solve(ocp, grid_size=30, print_level=0, tol=1e-8, init=init_constant)
65+
@testset verbose = true showtiming = true ":goddard :max_rf :all_constraints_types :init_constant" begin
66+
@test sol2.objective 1.0125 rtol=1e-2
67+
end
68+
69+
70+
# with initial guess from solution
71+
init_sol = OptimalControlInit(sol2)
72+
sol3 = solve(ocp, grid_size=30, print_level=0, tol=1e-8, init=init_sol)
73+
@testset verbose = true showtiming = true ":goddard :max_rf :all_constraints_types :init_sol" begin
74+
@test sol3.objective 1.0125 rtol=1e-2
75+
end

test/suite/simple_integrator.jl

Lines changed: 17 additions & 2 deletions
Original file line numberDiff line numberDiff line change
@@ -12,8 +12,23 @@ constraint!(ocp1, :initial, -1, :initial_constraint)
1212
constraint!(ocp1, :final, 0, :final_constraint)
1313
dynamics!(ocp1, (x, u) -> -x + u)
1414
objective!(ocp1, :lagrange, (x, u) -> u^2)
15-
init_constant = [-0.5, 0]
16-
sol1 = solve(ocp1, grid_size=100, print_level=0, tol=1e-12, init=init_constant)
15+
sol1 = solve(ocp1, grid_size=100, print_level=0, tol=1e-12)
1716
@testset verbose = true showtiming = true ":double_integrator :min_tf" begin
1817
@test sol1.objective 0.313 rtol=1e-2
1918
end
19+
20+
# with initial guess (using both vector and scalar syntax, no optimization variables)
21+
x_init = [-0.5]
22+
u_init = 0
23+
init_constant = OptimalControlInit(x_init, u_init)
24+
sol2 = solve(ocp1, grid_size=100, print_level=0, tol=1e-12, init=init_constant)
25+
@testset verbose = true showtiming = true ":double_integrator :min_tf :init_constant" begin
26+
@test sol2.objective 0.313 rtol=1e-2
27+
end
28+
29+
# with initial guess from solution
30+
init_sol = OptimalControlInit(sol2)
31+
sol3 = solve(ocp1, grid_size=100, print_level=0, tol=1e-12, init=init_sol)
32+
@testset verbose = true showtiming = true ":double_integrator :min_tf :init_sol" begin
33+
@test sol3.objective 0.313 rtol=1e-2
34+
end

0 commit comments

Comments
 (0)