Note

You can download this example as a Jupyter notebook or start it in interactive mode.

Redispatch Example with SciGRID network#

In this example, we compare a 2-stage market with an initial market clearing in two bidding zones with flow-based market coupling and a subsequent redispatch market (incl. curtailment) to an idealised nodal pricing scheme.

[1]:
import pypsa
import matplotlib.pyplot as plt
import cartopy.crs as ccrs
from pypsa.descriptors import get_switchable_as_dense as as_dense
[2]:
solver = "cbc"

Load example network#

[3]:
o = pypsa.examples.scigrid_de(from_master=True)
o.lines.s_max_pu = 0.7
o.lines.loc[["316", "527", "602"], "s_nom"] = 1715
o.set_snapshots([o.snapshots[12]])
WARNING:pypsa.io:Importing network from PyPSA version v0.17.1 while current version is v0.23.0. Read the release notes at https://pypsa.readthedocs.io/en/latest/release_notes.html to prepare your network for import.
INFO:pypsa.io:Imported network scigrid-de.nc has buses, generators, lines, loads, storage_units, transformers
[4]:
n = o.copy()  # for redispatch model
m = o.copy()  # for market model
[5]:
o.plot();
/home/docs/checkouts/readthedocs.org/user_builds/pypsa/conda/latest/lib/python3.11/site-packages/cartopy/mpl/style.py:76: UserWarning: facecolor will have no effect as it has been defined as "never".
  warnings.warn('facecolor will have no effect as it has been '
../_images/examples_scigrid-redispatch_6_1.png

Solve original nodal market model o#

First, let us solve a nodal market using the original model o:

[6]:
o.optimize(solver_name=solver)
WARNING:pypsa.components:The following transformers have zero r, which could break the linear load flow:
Index(['2', '5', '10', '12', '13', '15', '18', '20', '22', '24', '26', '30',
       '32', '37', '42', '46', '52', '56', '61', '68', '69', '74', '78', '86',
       '87', '94', '95', '96', '99', '100', '104', '105', '106', '107', '117',
       '120', '123', '124', '125', '128', '129', '138', '143', '156', '157',
       '159', '160', '165', '184', '191', '195', '201', '220', '231', '232',
       '233', '236', '247', '248', '250', '251', '252', '261', '263', '264',
       '267', '272', '279', '281', '282', '292', '303', '307', '308', '312',
       '315', '317', '322', '332', '334', '336', '338', '351', '353', '360',
       '362', '382', '384', '385', '391', '403', '404', '413', '421', '450',
       '458'],
      dtype='object', name='Transformer')
WARNING:pypsa.components:The following transformers have zero r, which could break the linear load flow:
Index(['2', '5', '10', '12', '13', '15', '18', '20', '22', '24', '26', '30',
       '32', '37', '42', '46', '52', '56', '61', '68', '69', '74', '78', '86',
       '87', '94', '95', '96', '99', '100', '104', '105', '106', '107', '117',
       '120', '123', '124', '125', '128', '129', '138', '143', '156', '157',
       '159', '160', '165', '184', '191', '195', '201', '220', '231', '232',
       '233', '236', '247', '248', '250', '251', '252', '261', '263', '264',
       '267', '272', '279', '281', '282', '292', '303', '307', '308', '312',
       '315', '317', '322', '332', '334', '336', '338', '351', '353', '360',
       '362', '382', '384', '385', '391', '403', '404', '413', '421', '450',
       '458'],
      dtype='object', name='Transformer')
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
INFO:linopy.model: Solve linear problem using Cbc solver
INFO:linopy.io: Writing time: 0.23s
INFO:linopy.constants: Optimization successful:
Status: ok
Termination condition: optimal
Solution: 2485 primals, 5957 duals
Objective: 3.01e+05
Solver model: not available
Solver message: Optimal - objective value 301209.38232509


Welcome to the CBC MILP Solver
Version: 2.10.10
Build Date: Apr 19 2023

command line - cbc -printingOptions all -import /tmp/linopy-problem-_16o4ioi.lp -solve -solu /tmp/linopy-solve-jjy71qft.sol (default strategy 1)
Option for printingOptions changed from normal to all
Presolve 619 (-5338) rows, 1084 (-1401) columns and 4185 (-7071) elements
Perturbing problem by 0.001% of 2290.1846 - largest nonzero change 0.00099399927 ( 0.0065461274%) - largest zero change 0.00099328721
0  Obj -9.7707246 Primal inf 1265555.5 (571)
87  Obj -9.7707246 Primal inf 836719.61 (538)
174  Obj -9.2332508 Primal inf 651250.73 (494)
247  Obj -8.8448867 Primal inf 978852.73 (487)
326  Obj -6.7711985 Primal inf 417327.98 (394)
402  Obj -4.7279448 Primal inf 950934.8 (415)
478  Obj 3999.92 Primal inf 898091.71 (286)
554  Obj 4034.688 Primal inf 745280.01 (250)
627  Obj 119709.84 Primal inf 25480.829 (97)
714  Obj 299100.97 Primal inf 510.73235 (26)
739  Obj 301211.97
Optimal - objective value 301209.38
After Postsolve, objective 301209.38, infeasibilities - dual 24.116221 (1), primal 6.0436268e-07 (1)
Presolved model was optimal, full model needs cleaning up
0  Obj 301209.38 Dual inf 0.24116211 (1)
End of values pass after 1 iterations
1  Obj 301209.38
Optimal - objective value 301209.38
Optimal objective 301209.3823 - 740 iterations time 0.122, Presolve 0.02
Total time (CPU seconds):       0.17   (Wallclock seconds):       0.12

[6]:
('ok', 'optimal')

Costs are 301 k€.

Build market model m with two bidding zones#

For this example, we split the German transmission network into two market zones at latitude 51 degrees.

You can build any other market zones by providing an alternative mapping from bus to zone.

[7]:
zones = (n.buses.y > 51).map(lambda x: "North" if x else "South")

Next, we assign this mapping to the market model m.

We re-assign the buses of all generators and loads, and remove all transmission lines within each bidding zone.

Here, we assume that the bidding zones are coupled through the transmission lines that connect them.

[8]:
for c in m.iterate_components(m.one_port_components):
    c.df.bus = c.df.bus.map(zones)

for c in m.iterate_components(m.branch_components):
    c.df.bus0 = c.df.bus0.map(zones)
    c.df.bus1 = c.df.bus1.map(zones)
    internal = c.df.bus0 == c.df.bus1
    m.mremove(c.name, c.df.loc[internal].index)

m.mremove("Bus", m.buses.index)
m.madd("Bus", ["North", "South"]);

Now, we can solve the coupled market with two bidding zones.

[9]:
m.optimize(solver_name=solver)
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
INFO:linopy.model: Solve linear problem using Cbc solver
INFO:linopy.io: Writing time: 0.19s
INFO:linopy.constants: Optimization successful:
Status: ok
Termination condition: optimal
Solution: 1561 primals, 3185 duals
Objective: 2.14e+05
Solver model: not available
Solver message: Optimal - objective value 213988.68595810


Welcome to the CBC MILP Solver
Version: 2.10.10
Build Date: Apr 19 2023

command line - cbc -printingOptions all -import /tmp/linopy-problem-ck53w22h.lp -solve -solu /tmp/linopy-solve-v67vcohp.sol (default strategy 1)
Option for printingOptions changed from normal to all
Presolve 40 (-3145) rows, 410 (-1151) columns and 487 (-4342) elements
Perturbing problem by 0.001% of 212.59539 - largest nonzero change 0.00017578427 ( 0.0036987348%) - largest zero change 0.00015445146
0  Obj 0 Primal inf 11285.222 (1)
48  Obj 184184.9 Primal inf 1700.1029 (24)
86  Obj 213988.73
Optimal - objective value 213988.69
After Postsolve, objective 213988.69, infeasibilities - dual 0 (0), primal 0 (0)
Optimal objective 213988.686 - 86 iterations time 0.012, Presolve 0.01
Total time (CPU seconds):       0.05   (Wallclock seconds):       0.03

[9]:
('ok', 'optimal')

Costs are 214 k€, which is much lower than the 301 k€ of the nodal market.

This is because network restrictions apart from the North/South division are not taken into account yet.

We can look at the market clearing prices of each zone:

[10]:
m.buses_t.marginal_price
[10]:
Bus North South
snapshot
2011-01-01 12:00:00 8.0 25.0

Build redispatch model n#

Next, based on the market outcome with two bidding zones m, we build a secondary redispatch market n that rectifies transmission constraints through curtailment and ramping up/down thermal generators.

First, we fix the dispatch of generators to the results from the market simulation. (For simplicity, this example disregards storage units.)

[11]:
p = m.generators_t.p / m.generators.p_nom
n.generators_t.p_min_pu = p
n.generators_t.p_max_pu = p

Then, we add generators bidding into redispatch market using the following assumptions:

  • All generators can reduce their dispatch to zero. This includes also curtailment of renewables.

  • All generators can increase their dispatch to their available/nominal capacity.

  • No changes to the marginal costs, i.e. reducing dispatch lowers costs.

With these settings, the 2-stage market should result in the same cost as the nodal market.

[12]:
g_up = n.generators.copy()
g_down = n.generators.copy()

g_up.index = g_up.index.map(lambda x: x + " ramp up")
g_down.index = g_down.index.map(lambda x: x + " ramp down")

up = (
    as_dense(m, "Generator", "p_max_pu") * m.generators.p_nom - m.generators_t.p
).clip(0) / m.generators.p_nom
down = -m.generators_t.p / m.generators.p_nom

up.columns = up.columns.map(lambda x: x + " ramp up")
down.columns = down.columns.map(lambda x: x + " ramp down")

n.madd("Generator", g_up.index, p_max_pu=up, **g_up.drop("p_max_pu", axis=1))

n.madd(
    "Generator",
    g_down.index,
    p_min_pu=down,
    p_max_pu=0,
    **g_down.drop(["p_max_pu", "p_min_pu"], axis=1)
);

Now, let’s solve the redispatch market:

[13]:
n.optimize(solver_name=solver)
WARNING:pypsa.components:The following transformers have zero r, which could break the linear load flow:
Index(['2', '5', '10', '12', '13', '15', '18', '20', '22', '24', '26', '30',
       '32', '37', '42', '46', '52', '56', '61', '68', '69', '74', '78', '86',
       '87', '94', '95', '96', '99', '100', '104', '105', '106', '107', '117',
       '120', '123', '124', '125', '128', '129', '138', '143', '156', '157',
       '159', '160', '165', '184', '191', '195', '201', '220', '231', '232',
       '233', '236', '247', '248', '250', '251', '252', '261', '263', '264',
       '267', '272', '279', '281', '282', '292', '303', '307', '308', '312',
       '315', '317', '322', '332', '334', '336', '338', '351', '353', '360',
       '362', '382', '384', '385', '391', '403', '404', '413', '421', '450',
       '458'],
      dtype='object', name='Transformer')
WARNING:pypsa.components:The following transformers have zero r, which could break the linear load flow:
Index(['2', '5', '10', '12', '13', '15', '18', '20', '22', '24', '26', '30',
       '32', '37', '42', '46', '52', '56', '61', '68', '69', '74', '78', '86',
       '87', '94', '95', '96', '99', '100', '104', '105', '106', '107', '117',
       '120', '123', '124', '125', '128', '129', '138', '143', '156', '157',
       '159', '160', '165', '184', '191', '195', '201', '220', '231', '232',
       '233', '236', '247', '248', '250', '251', '252', '261', '263', '264',
       '267', '272', '279', '281', '282', '292', '303', '307', '308', '312',
       '315', '317', '322', '332', '334', '336', '338', '351', '353', '360',
       '362', '382', '384', '385', '391', '403', '404', '413', '421', '450',
       '458'],
      dtype='object', name='Transformer')
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
INFO:linopy.model: Solve linear problem using Cbc solver
INFO:linopy.io: Writing time: 0.26s
INFO:linopy.constants: Optimization successful:
Status: ok
Termination condition: optimal
Solution: 5331 primals, 11649 duals
Objective: 3.01e+05
Solver model: not available
Solver message: Optimal - objective value 301209.38114435


Welcome to the CBC MILP Solver
Version: 2.10.10
Build Date: Apr 19 2023

command line - cbc -printingOptions all -import /tmp/linopy-problem-try2gvya.lp -solve -solu /tmp/linopy-solve-inay2oj2.sol (default strategy 1)
Option for printingOptions changed from normal to all
Presolve 624 (-11025) rows, 1312 (-4019) columns and 4432 (-15362) elements
Perturbing problem by 0.001% of 2290.1846 - largest nonzero change 0.00098855718 ( 0.0078811671%) - largest zero change 0.00098054313
0  Obj 194177.55 Primal inf 1293350.5 (574) Dual inf 7119.4271 (158)
87  Obj -9.7647787 Primal inf 877325.61 (543)
174  Obj -9.5234671 Primal inf 1640907.1 (537)
261  Obj -8.3005102 Primal inf 2420805.7 (518)
348  Obj -7.3783209 Primal inf 2147254 (477)
418  Obj -5.2529737 Primal inf 3244639.9 (438)
490  Obj 3998.7305 Primal inf 1175337.1 (376)
560  Obj 4000.8315 Primal inf 377613.76 (330)
642  Obj 4042.578 Primal inf 38545509 (362)
729  Obj 182923.2 Primal inf 11866.249 (80)
816  Obj 301201.3 Primal inf 0.96385324 (4)
819  Obj 301211.15
819  Obj 301209.38 Dual inf 1.4139822e-05 (3)
822  Obj 301209.38
Optimal - objective value 301209.38
After Postsolve, objective 301209.38, infeasibilities - dual 1473.5246 (101), primal 2.2530437e-05 (94)
Presolved model was optimal, full model needs cleaning up
0  Obj 301209.38 Primal inf 8.66223e-07 (5) Dual inf 5.0000001e+08 (106)
End of values pass after 108 iterations
108  Obj 301209.38
Optimal - objective value 301209.38
Optimal objective 301209.3811 - 930 iterations time 0.122, Presolve 0.03
Total time (CPU seconds):       0.21   (Wallclock seconds):       0.18

[13]:
('ok', 'optimal')

And, as expected, the costs are the same as for the nodal market: 301 k€.

Now, we can plot both the market results of the 2 bidding zone market and the redispatch results:

[14]:
fig, axs = plt.subplots(
    1, 3, figsize=(20, 10), subplot_kw={"projection": ccrs.AlbersEqualArea()}
)

market = (
    n.generators_t.p[m.generators.index]
    .T.squeeze()
    .groupby(n.generators.bus)
    .sum()
    .div(2e4)
)
n.plot(ax=axs[0], bus_sizes=market, title="2 bidding zones market simulation")

redispatch_up = (
    n.generators_t.p.filter(like="ramp up")
    .T.squeeze()
    .groupby(n.generators.bus)
    .sum()
    .div(2e4)
)
n.plot(
    ax=axs[1], bus_sizes=redispatch_up, bus_colors="blue", title="Redispatch: ramp up"
)

redispatch_down = (
    n.generators_t.p.filter(like="ramp down")
    .T.squeeze()
    .groupby(n.generators.bus)
    .sum()
    .div(-2e4)
)
n.plot(
    ax=axs[2],
    bus_sizes=redispatch_down,
    bus_colors="red",
    title="Redispatch: ramp down / curtail",
);
/home/docs/checkouts/readthedocs.org/user_builds/pypsa/conda/latest/lib/python3.11/site-packages/cartopy/mpl/style.py:76: UserWarning: facecolor will have no effect as it has been defined as "never".
  warnings.warn('facecolor will have no effect as it has been '
../_images/examples_scigrid-redispatch_30_1.png

We can also read out the final dispatch of each generator:

[15]:
grouper = n.generators.index.str.split(" ramp", expand=True).get_level_values(0)

n.generators_t.p.groupby(grouper, axis=1).sum().squeeze()
[15]:
1 Gas                     0.000000
1 Hard Coal               0.000000
1 Solar                  11.326192
1 Wind Onshore            1.754375
100_220kV Solar          14.913326
                           ...
98 Wind Onshore          71.451646
99_220kV Gas              0.000000
99_220kV Hard Coal        0.000000
99_220kV Solar            8.246606
99_220kV Wind Onshore     3.432939
Name: 2011-01-01 12:00:00, Length: 1423, dtype: float64

Changing bidding strategies in redispatch market#

We can also formulate other bidding strategies or compensation mechanisms for the redispatch market.

For example, that ramping up a generator is twice as expensive.

[16]:
n.generators.loc[n.generators.index.str.contains("ramp up"), "marginal_cost"] *= 2

Or that generators need to be compensated for curtailing them or ramping them down at 50% of their marginal cost.

[17]:
n.generators.loc[n.generators.index.str.contains("ramp down"), "marginal_cost"] *= -0.5

In this way, the outcome should be more expensive than the ideal nodal market:

[18]:
n.optimize(solver_name=solver)
WARNING:pypsa.components:The following transformers have zero r, which could break the linear load flow:
Index(['2', '5', '10', '12', '13', '15', '18', '20', '22', '24', '26', '30',
       '32', '37', '42', '46', '52', '56', '61', '68', '69', '74', '78', '86',
       '87', '94', '95', '96', '99', '100', '104', '105', '106', '107', '117',
       '120', '123', '124', '125', '128', '129', '138', '143', '156', '157',
       '159', '160', '165', '184', '191', '195', '201', '220', '231', '232',
       '233', '236', '247', '248', '250', '251', '252', '261', '263', '264',
       '267', '272', '279', '281', '282', '292', '303', '307', '308', '312',
       '315', '317', '322', '332', '334', '336', '338', '351', '353', '360',
       '362', '382', '384', '385', '391', '403', '404', '413', '421', '450',
       '458'],
      dtype='object', name='Transformer')
WARNING:pypsa.components:The following transformers have zero r, which could break the linear load flow:
Index(['2', '5', '10', '12', '13', '15', '18', '20', '22', '24', '26', '30',
       '32', '37', '42', '46', '52', '56', '61', '68', '69', '74', '78', '86',
       '87', '94', '95', '96', '99', '100', '104', '105', '106', '107', '117',
       '120', '123', '124', '125', '128', '129', '138', '143', '156', '157',
       '159', '160', '165', '184', '191', '195', '201', '220', '231', '232',
       '233', '236', '247', '248', '250', '251', '252', '261', '263', '264',
       '267', '272', '279', '281', '282', '292', '303', '307', '308', '312',
       '315', '317', '322', '332', '334', '336', '338', '351', '353', '360',
       '362', '382', '384', '385', '391', '403', '404', '413', '421', '450',
       '458'],
      dtype='object', name='Transformer')
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
<__array_function__ internals>:200: RuntimeWarning: invalid value encountered in cast
INFO:linopy.model: Solve linear problem using Cbc solver
INFO:linopy.io: Writing time: 0.28s
INFO:linopy.constants: Optimization successful:
Status: ok
Termination condition: optimal
Solution: 5331 primals, 11649 duals
Objective: 4.79e+05
Solver model: not available
Solver message: Optimal - objective value 479003.12190570


Welcome to the CBC MILP Solver
Version: 2.10.10
Build Date: Apr 19 2023

command line - cbc -printingOptions all -import /tmp/linopy-problem-1xxbdwqs.lp -solve -solu /tmp/linopy-solve-kif56wqi.sol (default strategy 1)
Option for printingOptions changed from normal to all
Presolve 624 (-11025) rows, 1314 (-4017) columns and 4435 (-15359) elements
Perturbing problem by 0.001% of 4580.3693 - largest nonzero change 0.00072687858 ( 0.0018687491%) - largest zero change 0.00072396771
0  Obj 223877.48 Primal inf 1293350.5 (574)
87  Obj 223877.76 Primal inf 1018408.7 (542)
174  Obj 223878.4 Primal inf 1372137.1 (529)
261  Obj 223879.83 Primal inf 1791131.4 (491)
340  Obj 223881.63 Primal inf 996329.82 (435)
407  Obj 223883.12 Primal inf 1172195 (381)
485  Obj 231268.34 Primal inf 2418876.4 (400)
556  Obj 231720.83 Primal inf 87863.979 (177)
636  Obj 232716.78 Primal inf 2722708.1 (299)
715  Obj 243559.7 Primal inf 14877.568 (59)
802  Obj 472596.31 Primal inf 51058.681 (61)
832  Obj 479005.19
832  Obj 479003.13 Dual inf 0.00029016411 (7)
839  Obj 479003.12
Optimal - objective value 479003.12
After Postsolve, objective 479003.12, infeasibilities - dual 2771.9398 (91), primal 1.9954696e-05 (85)
Presolved model was optimal, full model needs cleaning up
0  Obj 479003.12 Primal inf 7.2767925e-07 (4) Dual inf 4.0000003e+08 (95)
End of values pass after 96 iterations
96  Obj 479003.12
Optimal - objective value 479003.12
Optimal objective 479003.1219 - 935 iterations time 0.112, Presolve 0.03
Total time (CPU seconds):       0.22   (Wallclock seconds):       0.18

[18]:
('ok', 'optimal')

Costs are now 502 k€ compared to 301 k€.