Pricing callback
The pricing callback lets you define how to solve the subproblems of a Dantzig-Wolfe decomposition to generate a new entering column in the master program. This callback is useful when you know an efficient algorithm to solve the subproblems, i.e. an algorithm better than solving the subproblem with a MIP solver.
First, we load the packages and define aliases :
using Coluna, BlockDecomposition, JuMP, MathOptInterface, GLPK;
const BD = BlockDecomposition;
const MOI = MathOptInterface;
Let us see an example with the following generalized assignment problem :
M = 1:3;
J = 1:15;
c = [12.7 22.5 8.9 20.8 13.6 12.4 24.8 19.1 11.5 17.4 24.7 6.8 21.7 14.3 10.5; 19.1 24.8 24.4 23.6 16.1 20.6 15.0 9.5 7.9 11.3 22.6 8.0 21.5 14.7 23.2; 18.6 14.1 22.7 9.9 24.2 24.5 20.8 12.9 17.7 11.9 18.7 10.1 9.1 8.9 7.7; 13.1 16.2 16.8 16.7 9.0 16.9 17.9 12.1 17.5 22.0 19.9 14.6 18.2 19.6 24.2];
w = [61 70 57 82 51 74 98 64 86 80 69 79 60 76 78; 50 57 61 83 81 79 63 99 82 59 83 91 59 99 91;91 81 66 63 59 81 87 90 65 55 57 68 92 91 86; 62 79 73 60 75 66 68 99 69 60 56 100 67 68 54];
Q = [1020 1460 1530];
with the following Coluna configuration
coluna = optimizer_with_attributes(
Coluna.Optimizer,
"params" => Coluna.Params(
solver = Coluna.Algorithm.TreeSearchAlgorithm() # default BCP
),
"default_optimizer" => GLPK.Optimizer # GLPK for the master & the subproblems
);
for which the JuMP model takes the form:
model = BlockModel(coluna);
@axis(M_axis, M);
@variable(model, x[m in M_axis, j in J], Bin);
@constraint(model, cov[j in J], sum(x[m,j] for m in M_axis) == 1);
@objective(model, Min, sum(c[m,j]*x[m,j] for m in M_axis, j in J));
@dantzig_wolfe_decomposition(model, dwdec, M_axis);
where, as you can see, we omitted the knapsack constraints. These constraints are implicitly defined by the algorithm called in the pricing callback.
Let's use a knapsack algorithm defined by the following function to solve the knapsack subproblems:
function solve_knapsack(cost, weight, capacity)
sp_model = Model(GLPK.Optimizer)
items = 1:length(weight)
@variable(sp_model, x[i in items], Bin)
@constraint(sp_model, weight' * x <= capacity)
@objective(sp_model, Min, cost' * x)
optimize!(sp_model)
x_val = value.(x)
return filter(i -> x_val[i] ≈ 1, collect(items))
end
solve_knapsack (generic function with 1 method)
You can replace the content of the function with any algorithm that solves the knapsack problem (such as algorithms provided by the unregistered package Knapsacks).
The pricing callback is a function. It takes as argument cbdata
which is a data structure that allows the user to interact with Coluna within the pricing callback.
function my_pricing_callback(cbdata)
# Retrieve the index of the subproblem (it will be one of the values in M_axis)
cur_machine = BD.callback_spid(cbdata, model)
# Uncomment to see that the pricing callback is called.
# println("Pricing callback for machine $(cur_machine).")
# Retrieve reduced costs of subproblem variables
red_costs = [BD.callback_reduced_cost(cbdata, x[cur_machine, j]) for j in J]
# Run the knapsack algorithm
jobs_assigned_to_cur_machine = solve_knapsack(red_costs, w[cur_machine, :], Q[cur_machine])
# Create the solution (send only variables with non-zero values)
sol_vars = VariableRef[x[cur_machine, j] for j in jobs_assigned_to_cur_machine]
sol_vals = [1.0 for _ in jobs_assigned_to_cur_machine]
sol_cost = sum(red_costs[j] for j in jobs_assigned_to_cur_machine; init = 0.0)
# Submit the solution to the subproblem to Coluna
MOI.submit(model, BD.PricingSolution(cbdata), sol_cost, sol_vars, sol_vals)
# Submit the dual bound to the solution of the subproblem
# This bound is used to compute the contribution of the subproblem to the lagrangian
# bound in column generation.
MOI.submit(model, BD.PricingDualBound(cbdata), sol_cost) # optimal solution
return
end
my_pricing_callback (generic function with 1 method)
The pricing callback is provided to Coluna using the keyword solver
in the method specify!
.
subproblems = BD.getsubproblems(dwdec);
BD.specify!.(subproblems, lower_multiplicity = 0, solver = my_pricing_callback);
You can then optimize :
optimize!(model);
Coluna
Version 0.8.2 | https://github.com/atoptima/Coluna.jl
***************************************************************************************
**** B&B tree root node
**** Local DB = -Inf, global bounds: [ -Inf , Inf ], time = 0.00 sec.
***************************************************************************************
<st= 1> <it= 1> <et= 0.38> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-199514.2000> <mlp=100000.0000> <PB=Inf>
<st= 1> <it= 2> <et= 0.38> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-249337.6000> <mlp=50148.8000> <PB=Inf>
<st= 1> <it= 3> <et= 0.39> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-249373.5000> <mlp=50145.2000> <PB=Inf>
<st= 1> <it= 4> <et= 0.39> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-180037.2000> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 5> <et= 0.39> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-179876.9500> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 6> <et= 0.39> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-179800.0500> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 7> <et= 0.39> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-179623.2500> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 8> <et= 0.39> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-172484.4875> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 9> <et= 0.40> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-179618.9000> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 10> <et= 0.40> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-164410.0125> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 11> <et= 0.40> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-147610.4769> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 12> <et= 0.40> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-134418.0173> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 13> <et= 0.40> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-107799.0925> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 14> <et= 0.41> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB=-92017.9579> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 15> <et= 0.41> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= -14.7849> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 16> <et= 0.41> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 157.8132> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 17> <et= 0.41> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= -14.7699> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 18> <et= 0.41> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= -315.3353> <mlp= 216.6500> <PB=Inf>
<st= 1> <it= 19> <et= 0.41> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= -328.0929> <mlp= 212.3356> <PB=Inf>
<st= 1> <it= 20> <et= 0.42> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 90.9151> <mlp= 199.2035> <PB=Inf>
<st= 1> <it= 21> <et= 0.42> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 99.4029> <mlp= 194.8235> <PB=Inf>
<st= 1> <it= 22> <et= 0.42> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 104.8426> <mlp= 191.6053> <PB=Inf>
<st= 1> <it= 23> <et= 0.42> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 118.2290> <mlp= 189.6338> <PB=Inf>
<st= 1> <it= 24> <et= 0.42> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 77.5000> <mlp= 183.4000> <PB=Inf>
<st= 1> <it= 25> <et= 0.43> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 112.2081> <mlp= 182.0297> <PB=Inf>
<st= 1> <it= 26> <et= 0.43> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 122.0606> <mlp= 181.0152> <PB=Inf>
<st= 1> <it= 27> <et= 0.43> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 104.9603> <mlp= 179.9977> <PB=Inf>
<st= 1> <it= 28> <et= 0.43> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 110.5347> <mlp= 179.0653> <PB=Inf>
<st= 1> <it= 29> <et= 0.43> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 149.6000> <mlp= 178.0745> <PB=Inf>
<st= 1> <it= 30> <et= 0.44> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 125.2250> <mlp= 176.8750> <PB=Inf>
<st= 1> <it= 31> <et= 0.44> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 154.2188> <mlp= 176.4187> <PB=Inf>
<st= 1> <it= 32> <et= 0.44> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 138.0000> <mlp= 171.8500> <PB=Inf>
<st= 1> <it= 33> <et= 0.44> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 158.2400> <mlp= 171.8500> <PB=Inf>
<st= 1> <it= 34> <et= 0.45> <mst= 0.00> <sp= 0.00> <cols= 2> <al= 0.00> <DB= 159.4500> <mlp= 171.8500> <PB=Inf>
<st= 1> <it= 35> <et= 0.45> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 156.3000> <mlp= 171.8500> <PB=Inf>
<st= 1> <it= 36> <et= 0.45> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 162.0000> <mlp= 171.8500> <PB=Inf>
<st= 1> <it= 37> <et= 0.45> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 161.5500> <mlp= 171.7750> <PB=Inf>
<st= 1> <it= 38> <et= 0.45> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 163.6000> <mlp= 171.7750> <PB=Inf>
<st= 1> <it= 39> <et= 0.46> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 158.6000> <mlp= 169.7000> <PB=Inf>
[ Info: Improving primal solution with value 167.70000000000002 is found during column generation
<st= 1> <it= 40> <et= 0.46> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 164.7000> <mlp= 167.7000> <PB=167.7000>
<st= 1> <it= 41> <et= 0.46> <mst= 0.00> <sp= 0.00> <cols= 3> <al= 0.00> <DB= 164.7000> <mlp= 167.7000> <PB=167.7000>
<st= 1> <it= 42> <et= 0.46> <mst= 0.00> <sp= 0.00> <cols= 1> <al= 0.00> <DB= 166.5000> <mlp= 167.7000> <PB=167.7000>
<st= 1> <it= 43> <et= 0.46> <mst= 0.00> <sp= 0.00> <cols= 1> <al= 0.00> <DB= 166.5000> <mlp= 167.7000> <PB=167.7000>
[ Info: Improving primal solution with value 167.09999999999997 is found during column generation
<st= 1> <it= 44> <et= 0.47> <mst= 0.00> <sp= 0.00> <cols= 2> <al= 0.00> <DB= 166.1000> <mlp= 167.1000> <PB=167.1000>
<st= 1> <it= 45> <et= 0.47> <mst= 0.00> <sp= 0.00> <cols= 2> <al= 0.00> <DB= 166.1000> <mlp= 167.1000> <PB=167.1000>
<st= 1> <it= 46> <et= 0.47> <mst= 0.00> <sp= 0.00> <cols= 2> <al= 0.00> <DB= 165.7000> <mlp= 167.1000> <PB=167.1000>
<st= 1> <it= 47> <et= 0.47> <mst= 0.00> <sp= 0.00> <cols= 1> <al= 0.00> <DB= 166.5000> <mlp= 167.1000> <PB=167.1000>
[ Info: Improving primal solution with value 166.5 is found during column generation
<st= 1> <it= 48> <et= 0.47> <mst= 0.00> <sp= 0.00> <cols= 0> <al= 0.00> <DB= 166.5000> <mlp= 166.5000> <PB=166.5000>
───────────────────────────────────────────────────────────────────────────────
Time Allocations
─────────────────────── ────────────────────────
Tot / % measured: 1.05s / 45.1% 140MiB / 67.8%
Section ncalls time %tot avg alloc %tot avg
───────────────────────────────────────────────────────────────────────────────
Coluna 1 473ms 100.0% 473ms 95.1MiB 100.0% 95.1MiB
Pricing callback 144 451ms 95.4% 3.13ms 77.6MiB 81.6% 552KiB
SolveLpForm 48 12.2ms 2.6% 254μs 5.85MiB 6.2% 125KiB
───────────────────────────────────────────────────────────────────────────────
[ Info: Terminated
[ Info: Primal bound: 166.5
[ Info: Dual bound: 166.50000000000006
and retrieve the information you need as usual :
objective_value(model)
166.5
This page was generated using Literate.jl.