-
Notifications
You must be signed in to change notification settings - Fork 12
/
nlpmodels.jl
executable file
·542 lines (468 loc) · 16.3 KB
/
nlpmodels.jl
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
#!/usr/bin/env julia
###### AC-OPF using ADNLPModels ######
#
# implementation reference: https://juliasmoothoptimizers.github.io/ADNLPModels.jl/stable/tutorial/
# other AD libraries can be considered: https://juliasmoothoptimizers.github.io/ADNLPModels.jl/stable/
#
# The original rosetta-opf implementation can be found here,
# https://github.com/lanl-ansi/rosetta-opf/blob/bdc924c23694d09f6393575f6d732176aef7605d/nlpmodels.jl
#
# This implementation is an adaption of the original rosetta-opf implementation
# that uses ConcreteStructs based on the SciMLBenchmarks from 01/16/2024,
# https://github.com/SciML/SciMLBenchmarks.jl/blob/master/benchmarks/OptimizationFrameworks/optimal_powerflow.jmd
#
import PowerModels
import ADNLPModels
import ConcreteStructs
import NLPModelsIpopt
ConcreteStructs.@concrete struct DataRepresentation
data
ref
var_lookup
var_init
var_lb
var_ub
ref_gen_idxs
lookup_pg
lookup_qg
lookup_va
lookup_vm
lookup_lij
lookup_p_lij
lookup_q_lij
cost_arrs
f_bus
t_bus
ref_bus_idxs
ref_buses_idxs
ref_bus_gens
ref_bus_arcs
ref_branch_idxs
ref_arcs_from
ref_arcs_to
p_idxmap
q_idxmap
bus_pd
bus_qd
bus_gs
bus_bs
br_g
br_b
br_tr
br_ti
br_ttm
br_g_fr
br_b_fr
br_g_to
br_b_to
end
function load_and_setup_data(file_name)
data = PowerModels.parse_file(file_name)
PowerModels.standardize_cost_terms!(data, order=2)
PowerModels.calc_thermal_limits!(data)
ref = PowerModels.build_ref(data)[:it][:pm][:nw][0]
# Some data munging to type-stable forms
var_lookup = Dict{String,Int}()
var_init = Float64[]
var_lb = Float64[]
var_ub = Float64[]
var_idx = 1
for (i,bus) in ref[:bus]
push!(var_init, 0.0) #va
push!(var_lb, -Inf)
push!(var_ub, Inf)
var_lookup["va_$(i)"] = var_idx
var_idx += 1
push!(var_init, 1.0) #vm
push!(var_lb, bus["vmin"])
push!(var_ub, bus["vmax"])
var_lookup["vm_$(i)"] = var_idx
var_idx += 1
end
for (i,gen) in ref[:gen]
push!(var_init, 0.0) #pg
push!(var_lb, gen["pmin"])
push!(var_ub, gen["pmax"])
var_lookup["pg_$(i)"] = var_idx
var_idx += 1
push!(var_init, 0.0) #qg
push!(var_lb, gen["qmin"])
push!(var_ub, gen["qmax"])
var_lookup["qg_$(i)"] = var_idx
var_idx += 1
end
for (l,i,j) in ref[:arcs]
branch = ref[:branch][l]
push!(var_init, 0.0) #p
push!(var_lb, -branch["rate_a"])
push!(var_ub, branch["rate_a"])
var_lookup["p_$(l)_$(i)_$(j)"] = var_idx
var_idx += 1
push!(var_init, 0.0) #q
push!(var_lb, -branch["rate_a"])
push!(var_ub, branch["rate_a"])
var_lookup["q_$(l)_$(i)_$(j)"] = var_idx
var_idx += 1
end
@assert var_idx == length(var_init)+1
ref_gen_idxs = [i for i in keys(ref[:gen])]
lookup_pg = Dict{Int,Int}()
lookup_qg = Dict{Int,Int}()
lookup_va = Dict{Int,Int}()
lookup_vm = Dict{Int,Int}()
lookup_lij = Tuple{Int,Int,Int}[]
lookup_p_lij = Int[]
lookup_q_lij = Int[]
cost_arrs = Dict{Int,Vector{Float64}}()
for (i,gen) in ref[:gen]
lookup_pg[i] = var_lookup["pg_$(i)"]
lookup_qg[i] = var_lookup["qg_$(i)"]
cost_arrs[i] = gen["cost"]
end
for (i,bus) in ref[:bus]
lookup_va[i] = var_lookup["va_$(i)"]
lookup_vm[i] = var_lookup["vm_$(i)"]
end
for (l,i,j) in ref[:arcs]
push!(lookup_lij, (l,i,j))
push!(lookup_p_lij,var_lookup["p_$(l)_$(i)_$(j)"])
push!(lookup_q_lij,var_lookup["q_$(l)_$(i)_$(j)"])
end
f_bus = Dict{Int,Int}()
t_bus = Dict{Int,Int}()
for (l,branch) in ref[:branch]
f_bus[l] = branch["f_bus"]
t_bus[l] = branch["t_bus"]
end
ref_bus_idxs = [i for i in keys(ref[:bus])]
ref_buses_idxs = [i for i in keys(ref[:ref_buses])]
ref_bus_gens = ref[:bus_gens]
ref_bus_arcs = ref[:bus_arcs]
ref_branch_idxs = [i for i in keys(ref[:branch])]
ref_arcs_from = ref[:arcs_from]
ref_arcs_to = ref[:arcs_to]
p_idxmap = Dict(lookup_lij[i] => lookup_p_lij[i] for i in 1:length(lookup_lij))
q_idxmap = Dict(lookup_lij[i] => lookup_q_lij[i] for i in 1:length(lookup_lij))
bus_pd = Dict(i => 0.0 for (i,bus) in ref[:bus])
bus_qd = Dict(i => 0.0 for (i,bus) in ref[:bus])
bus_gs = Dict(i => 0.0 for (i,bus) in ref[:bus])
bus_bs = Dict(i => 0.0 for (i,bus) in ref[:bus])
for (i,bus) in ref[:bus]
if length(ref[:bus_loads][i]) > 0
bus_pd[i] = sum(ref[:load][l]["pd"] for l in ref[:bus_loads][i])
bus_qd[i] = sum(ref[:load][l]["qd"] for l in ref[:bus_loads][i])
end
if length(ref[:bus_shunts][i]) > 0
bus_gs[i] = sum(ref[:shunt][s]["gs"] for s in ref[:bus_shunts][i])
bus_bs[i] = sum(ref[:shunt][s]["bs"] for s in ref[:bus_shunts][i])
end
end
br_g = Dict(i => 0.0 for (i,branch) in ref[:branch])
br_b = Dict(i => 0.0 for (i,branch) in ref[:branch])
br_tr = Dict(i => 0.0 for (i,branch) in ref[:branch])
br_ti = Dict(i => 0.0 for (i,branch) in ref[:branch])
br_ttm = Dict(i => 0.0 for (i,branch) in ref[:branch])
br_g_fr = Dict(i => 0.0 for (i,branch) in ref[:branch])
br_b_fr = Dict(i => 0.0 for (i,branch) in ref[:branch])
br_g_to = Dict(i => 0.0 for (i,branch) in ref[:branch])
br_b_to = Dict(i => 0.0 for (i,branch) in ref[:branch])
for (i,branch) in ref[:branch]
g, b = PowerModels.calc_branch_y(branch)
tr, ti = PowerModels.calc_branch_t(branch)
br_g[i] = g
br_b[i] = b
br_tr[i] = tr
br_ti[i] = ti
br_ttm[i] = tr^2 + ti^2
br_g_fr[i] = branch["g_fr"]
br_b_fr[i] = branch["b_fr"]
br_g_to[i] = branch["g_to"]
br_b_to[i] = branch["b_to"]
end
return DataRepresentation(
data,
ref,
var_lookup,
var_init,
var_lb,
var_ub,
ref_gen_idxs,
lookup_pg,
lookup_qg,
lookup_va,
lookup_vm,
lookup_lij,
lookup_p_lij,
lookup_q_lij,
cost_arrs,
f_bus,
t_bus,
ref_bus_idxs,
ref_buses_idxs,
ref_bus_gens,
ref_bus_arcs,
ref_branch_idxs,
ref_arcs_from,
ref_arcs_to,
p_idxmap,
q_idxmap,
bus_pd,
bus_qd,
bus_gs,
bus_bs,
br_g,
br_b,
br_tr,
br_ti,
br_ttm,
br_g_fr,
br_b_fr,
br_g_to,
br_b_to)
end
function build_opf_optimization_prob(dataset)
(;data,
ref,
var_lookup,
var_init,
var_lb,
var_ub,
ref_gen_idxs,
lookup_pg,
lookup_qg,
lookup_va,
lookup_vm,
lookup_lij,
lookup_p_lij,
lookup_q_lij,
cost_arrs,
f_bus,
t_bus,
ref_bus_idxs,
ref_buses_idxs,
ref_bus_gens,
ref_bus_arcs,
ref_branch_idxs,
ref_arcs_from,
ref_arcs_to,
p_idxmap,
q_idxmap,
bus_pd,
bus_qd,
bus_gs,
bus_bs,
br_g,
br_b,
br_tr,
br_ti,
br_ttm,
br_g_fr,
br_b_fr,
br_g_to,
br_b_to) = dataset
function opf_objective(x)
cost = 0.0
for i in ref_gen_idxs
pg = x[lookup_pg[i]]
_cost_arr = cost_arrs[i]
cost += _cost_arr[1]*pg^2 + _cost_arr[2]*pg + _cost_arr[3]
end
return cost
end
function opf_constraints!(ret, x)
offsetidx = 0
# va_con
for (reti,i) in enumerate(ref_buses_idxs)
ret[reti + offsetidx] = x[lookup_va[i]]
end
offsetidx += length(ref_buses_idxs)
# @constraint(model,
# sum(p[a] for a in ref[:bus_arcs][i]) ==
# sum(pg[g] for g in ref_bus_gens[i]) -
# sum(load["pd"] for load in bus_loads) -
# sum(shunt["gs"] for shunt in bus_shunts)*x[lookup_vm[i]]^2
# )
for (reti,i) in enumerate(ref_bus_idxs)
ret[reti + offsetidx] = sum(x[lookup_pg[j]] for j in ref_bus_gens[i]; init=0.0) -
bus_pd[i] -
bus_gs[i]*x[lookup_vm[i]]^2 -
sum(x[p_idxmap[a]] for a in ref_bus_arcs[i])
end
offsetidx += length(ref_bus_idxs)
# @constraint(model,
# sum(q[a] for a in ref[:bus_arcs][i]) ==
# sum(x[lookup_qg[g]] for g in ref_bus_gens[i]) -
# sum(load["qd"] for load in bus_loads) +
# sum(shunt["bs"] for shunt in bus_shunts)*x[lookup_vm[i]]^2
# )
for (reti,i) in enumerate(ref_bus_idxs)
ret[reti + offsetidx] = sum(x[lookup_qg[j]] for j in ref_bus_gens[i]; init=0.0) -
bus_qd[i] +
bus_bs[i]*x[lookup_vm[i]]^2 -
sum(x[q_idxmap[a]] for a in ref_bus_arcs[i])
end
offsetidx += length(ref_bus_idxs)
# @NLconstraint(model, p_fr == (g+g_fr)/ttm*vm_fr^2 + (-g*tr+b*ti)/ttm*(vm_fr*vm_to*cos(va_fr-va_to)) + (-b*tr-g*ti)/ttm*(vm_fr*vm_to*sin(va_fr-va_to)) )
# power_flow_p_from_con =
for (reti,(l,i,j)) in enumerate(ref_arcs_from)
ret[reti + offsetidx] = (br_g[l]+br_g_fr[l])/br_ttm[l]*x[lookup_vm[f_bus[l]]]^2 +
(-br_g[l]*br_tr[l]+br_b[l]*br_ti[l])/br_ttm[l]*(x[lookup_vm[f_bus[l]]]*x[lookup_vm[t_bus[l]]]*cos(x[lookup_va[f_bus[l]]]-x[lookup_va[t_bus[l]]])) +
(-br_b[l]*br_tr[l]-br_g[l]*br_ti[l])/br_ttm[l]*(x[lookup_vm[f_bus[l]]]*x[lookup_vm[t_bus[l]]]*sin(x[lookup_va[f_bus[l]]]-x[lookup_va[t_bus[l]]])) -
x[p_idxmap[(l,i,j)]]
end
offsetidx += length(ref_arcs_from)
# @NLconstraint(model, p_to == (g+g_to)*vm_to^2 + (-g*tr-b*ti)/ttm*(vm_to*vm_fr*cos(va_to-va_fr)) + (-b*tr+g*ti)/ttm*(vm_to*vm_fr*sin(va_to-va_fr)) )
# power_flow_p_to_con
for (reti,(l,i,j)) in enumerate(ref_arcs_to)
ret[reti + offsetidx] = (br_g[l]+br_g_to[l])*x[lookup_vm[t_bus[l]]]^2 +
(-br_g[l]*br_tr[l]-br_b[l]*br_ti[l])/br_ttm[l]*(x[lookup_vm[t_bus[l]]]*x[lookup_vm[f_bus[l]]]*cos(x[lookup_va[t_bus[l]]]-x[lookup_va[f_bus[l]]])) +
(-br_b[l]*br_tr[l]+br_g[l]*br_ti[l])/br_ttm[l]*(x[lookup_vm[t_bus[l]]]*x[lookup_vm[f_bus[l]]]*sin(x[lookup_va[t_bus[l]]]-x[lookup_va[f_bus[l]]])) -
x[p_idxmap[(l,i,j)]]
end
offsetidx += length(ref_arcs_to)
# @NLconstraint(model, q_fr == -(b+b_fr)/ttm*vm_fr^2 - (-b*tr-g*ti)/ttm*(vm_fr*vm_to*cos(va_fr-va_to)) + (-g*tr+b*ti)/ttm*(vm_fr*vm_to*sin(va_fr-va_to)) )
# power_flow_q_from_con
for (reti,(l,i,j)) in enumerate(ref_arcs_from)
ret[reti + offsetidx] = -(br_b[l]+br_b_fr[l])/br_ttm[l]*x[lookup_vm[f_bus[l]]]^2 -
(-br_b[l]*br_tr[l]-br_g[l]*br_ti[l])/br_ttm[l]*(x[lookup_vm[f_bus[l]]]*x[lookup_vm[t_bus[l]]]*cos(x[lookup_va[f_bus[l]]]-x[lookup_va[t_bus[l]]])) +
(-br_g[l]*br_tr[l]+br_b[l]*br_ti[l])/br_ttm[l]*(x[lookup_vm[f_bus[l]]]*x[lookup_vm[t_bus[l]]]*sin(x[lookup_va[f_bus[l]]]-x[lookup_va[t_bus[l]]])) -
x[q_idxmap[(l,i,j)]]
end
offsetidx += length(ref_arcs_from)
# @NLconstraint(model, q_to == -(b+b_to)*vm_to^2 - (-b*tr+g*ti)/ttm*(vm_to*vm_fr*cos(va_to-va_fr)) + (-g*tr-b*ti)/ttm*(vm_to*vm_fr*sin(va_to-va_fr)) )
# power_flow_q_to_con
for (reti,(l,i,j)) in enumerate(ref_arcs_to)
ret[reti + offsetidx] = -(br_b[l]+br_b_to[l])*x[lookup_vm[t_bus[l]]]^2 -
(-br_b[l]*br_tr[l]+br_g[l]*br_ti[l])/br_ttm[l]*(x[lookup_vm[t_bus[l]]]*x[lookup_vm[f_bus[l]]]*cos(x[lookup_va[t_bus[l]]]-x[lookup_va[f_bus[l]]])) +
(-br_g[l]*br_tr[l]-br_b[l]*br_ti[l])/br_ttm[l]*(x[lookup_vm[t_bus[l]]]*x[lookup_vm[f_bus[l]]]*sin(x[lookup_va[t_bus[l]]]-x[lookup_va[f_bus[l]]])) -
x[q_idxmap[(l,i,j)]]
end
offsetidx += length(ref_arcs_to)
# @constraint(model, va_fr - va_to <= branch["angmax"])
# @constraint(model, va_fr - va_to >= branch["angmin"])
# power_flow_vad_con
for (reti,(l,i,j)) in enumerate(ref_arcs_from)
ret[reti + offsetidx] = x[lookup_va[f_bus[l]]] - x[lookup_va[t_bus[l]]]
end
offsetidx += length(ref_arcs_from)
# @constraint(model, p_fr^2 + q_fr^2 <= branch["rate_a"]^2)
# power_flow_mva_from_con
for (reti,(l,i,j)) in enumerate(ref_arcs_from)
ret[reti + offsetidx] = x[p_idxmap[(l,i,j)]]^2 + x[q_idxmap[(l,i,j)]]^2
end
offsetidx += length(ref_arcs_from)
# @constraint(model, p_to^2 + q_to^2 <= branch["rate_a"]^2)
# power_flow_mva_to_con
for (reti,(l,i,j)) in enumerate(ref_arcs_to)
ret[reti + offsetidx] = x[p_idxmap[(l,i,j)]]^2 + x[q_idxmap[(l,i,j)]]^2
end
offsetidx += length(ref_arcs_to)
@assert offsetidx == length(ret)
return ret
end
con_lbs = Float64[]
con_ubs = Float64[]
#@constraint(model, va[i] == 0)
for (i,bus) in ref[:ref_buses]
push!(con_lbs, 0.0)
push!(con_ubs, 0.0)
end
#power_balance_p_con
for (i,bus) in ref[:bus]
push!(con_lbs, 0.0)
push!(con_ubs, 0.0)
end
#power_balance_q_con
for (i,bus) in ref[:bus]
push!(con_lbs, 0.0)
push!(con_ubs, 0.0)
end
#power_flow_p_from_con
for (l,i,j) in ref[:arcs_from]
push!(con_lbs, 0.0)
push!(con_ubs, 0.0)
end
#power_flow_p_to_con
for (l,i,j) in ref[:arcs_to]
push!(con_lbs, 0.0)
push!(con_ubs, 0.0)
end
#power_flow_q_from_con
for (l,i,j) in ref[:arcs_from]
push!(con_lbs, 0.0)
push!(con_ubs, 0.0)
end
#power_flow_q_to_con
for (l,i,j) in ref[:arcs_to]
push!(con_lbs, 0.0)
push!(con_ubs, 0.0)
end
#power_flow_vad_con
for (l,i,j) in ref[:arcs_from]
branch = ref[:branch][l]
push!(con_lbs, branch["angmin"])
push!(con_ubs, branch["angmax"])
end
#power_flow_mva_from_con
for (l,i,j) in ref[:arcs_from]
branch = ref[:branch][l]
push!(con_lbs, -Inf)
push!(con_ubs, branch["rate_a"]^2)
end
#power_flow_mva_to_con
for (l,i,j) in ref[:arcs_to]
branch = ref[:branch][l]
push!(con_lbs, -Inf)
push!(con_ubs, branch["rate_a"]^2)
end
nlp = ADNLPModels.ADNLPModel!(opf_objective, var_init, var_lb, var_ub, opf_constraints!, con_lbs, con_ubs, backend = :optimized)
return (
nlp = nlp,
con_lbs = con_lbs,
con_ubs = con_ubs,
)
end
function solve_opf(file_name)
start = time()
dataset = load_and_setup_data(file_name);
data_load_time = time() - start
model = build_opf_optimization_prob(dataset)
nlp = model.nlp
model_build_time = time() - start - data_load_time
output = NLPModelsIpopt.ipopt(nlp)
cost = output.objective
feasible = (output.primal_feas <= 1e-6)
total_time = time() - start
solve_time = total_time - model_build_time - data_load_time
model_variables = length(dataset.var_init)
model_constraints = length(model.con_lbs)
println("variables: $(model_variables), $(length(dataset.var_lb)), $(length(dataset.var_ub))")
println("constraints: $(model_constraints), $(length(model.con_lbs)), $(length(model.con_ubs))")
println("")
println("\033[1mSummary\033[0m")
println(" case........: $(file_name)")
println(" variables...: $(model_variables)")
println(" constraints.: $(model_constraints)")
println(" feasible....: $(feasible)")
println(" cost........: $(round(Int, cost))")
println(" total time..: $(total_time)")
println(" data time.: $(data_load_time)")
println(" build time: $(model_build_time)")
println(" solve time: $(solve_time)")
println("")
return Dict(
"case" => file_name,
"variables" => model_variables,
"constraints" => model_constraints,
"feasible" => feasible,
"cost" => cost,
"time_total" => total_time,
"time_data" => data_load_time,
"time_build" => model_build_time,
"time_solve" => solve_time,
"solution" => Dict(k => output.solution[v] for (k, v) in dataset.var_lookup),
)
end
if isinteractive() == false
solve_opf("$(@__DIR__)/data/opf_warmup.m")
end