API for MaximumWeightTwoStageSpanningTree.jl

Public

MaximumWeightTwoStageSpanningTree.bendersMethod
function two_stage_spanning_tree_benders(inst::TwoStageSpanningTreeInstance;MILP_solver=GLPK.Optimizer, tol=0.000001)

computes an optimal solution (with tolerance tol) of the two stage spanning tree problem using a Benders approach.
source
MaximumWeightTwoStageSpanningTree.build_solve_and_encode_instance_as_maximum_weight_spanning_tree_single_scenario_layerMethod
function build_solve_and_encode_instance_as_maximum_weight_spanning_tree_single_scenario_layer(;
    grid_size=3,
    seed=0,
    nb_scenarios=10, 
    first_max=10, 
    second_max=20, 
    solver=lagrangian_heuristic_solver, 
    load_and_save=true
)

Builds a two stage spanning tree instance with a square grid graph of width grid_size, seed for the random number generator, nb_scenarios for the second stage, first_max and second_max as first and second stage maximum weight.

Solves it with solver

Encodes it for a pipeline with a maxium weight spanning tree for a two stage instance with single scenario second stage layer.

source
MaximumWeightTwoStageSpanningTree.column_generationMethod
function column_generation_two_stage_spanning_tree(inst::TwoStageSpanningTreeInstance;MILP_solver=GLPK.Optimizer,tol=0.00001)

Solves the dual of the linear relaxation of the two stage spanning tree MILP using a constraint generation.

returns objective_value, duals

source
MaximumWeightTwoStageSpanningTree.cut_generationMethod
function two_stage_spanning_tree_cut_generation(inst::TwoStageSpanningTreeInstance;MILP_solver=GLPK.Optimizer, separate_constraint_function=separate_forest_polytope_constraint_vertex_set_using_min_cut_MILP_formulation!, tol=0.000001, silent=false)

Solve the TwoStageSpanningTree instance inst using a cut generation approach.

separate_constraint_function indicates which method is used to separate the constraint on subsets of vertices. One option is available

  • separate_forest_polytope_constraint_vertex_set_using_min_cut_MILP_formulation!
source
MaximumWeightTwoStageSpanningTree.lagrangian_relaxationMethod
function lagrangian_relaxation(inst::TwoStageSpanningTreeInstance; nb_epochs=100)

runs
 - a subgradient algorithm during nb_epochs to solve the Lagrangian relaxation of inst (with `θ` initialized to 0)
 - a lagrangian heuristic on the resulting solution

returns lb , ub, forest, θ, training_losses with
- lb: value of the Lagrangian relaxation
- ub: value of the solution computed by the lagrangian heuristic
- forest: solution computed by the lagrangian heuristic
- θ: second stage problem value
- training_losses: vector with the value of the training losses along the subgradient algorithm
source
MaximumWeightTwoStageSpanningTree.test_models_on_data_setMethod

function testmodelsondataset( ; models, # Model name dataset, # Dataset obtained with build_or_load_spanning_tree_CO_layer_datasets() )

returns an array, each instance of which contains a tuple (inst, lb, UBs) where

  • inst is an instance
  • lb is the lower bound (lagrangian or optimal solution)
  • Ubs is a dictionnary with the ub (lagrangian or exaxt), and the performance of every model on the instance
source
MaximumWeightTwoStageSpanningTree.train_bbl_model_and_add_to_dict!Method

function trainbblmodelandaddtodict!( models # Dictionnary to which the result will be added ; traindataname, # Data name (to save on disk) traindata, # Data obtained with `buildorloadspanningtreeCOlayerdatasets()` nbDIRECTiterations=1000, # Nb iterations of the black box algorithm (DIRECT) perturbed=false, # Activate Gaussian perturbation nbperturbations=20, # Nb pertubation scenarios perturbationintensity=0.1, # Strength of the pertubation force_recompute=false # learn even if model available on disk )

source
MaximumWeightTwoStageSpanningTree.train_save_or_load_FYL_model!Method

function trainsaveorloadFYLmodel!( ; nbsamples, # nbsamples used in FYL perturbation traindataname, # name of the training set (used for saving or loading) traindata, # dataset of instance obtained with build_or_load_spanning_tree_CO_layer_datasets() nbepochs # nbepochs in the stochastic gradient descent )

source

Private

MaximumWeightTwoStageSpanningTree.BlackBoxLossType
struct BlackBoxLoss{A<:AbstractArray,I,P} <: Function
    nb_features::Int
    nb_samples::Int
    training_set::Vector{Tuple{A,I,Float64}}
    cost_pipeline::P     # Cost ∘ Decoder ∘ Maximizer: Outputs a solution
end

Encodes all the information for learning a problem by experience.

source
MaximumWeightTwoStageSpanningTree.BlackBoxLossMethod
function BlackBoxLoss(training_data,maximizer,decoder,cost_function;scaling_function=x->1.0)

Constructor for a BlackBoxLoss

Pipeline: x –GLM–> θ –Maximizer–> y –Decoder–> z

  • training_data: Vector(x,inst) where inst is an instance of the problem and x its features encoding
  • maximizer(θ,inst)
  • decoder(y,inst)
  • cost_function(z,inst)
  • scaling_function: order of magnitude of cost_function(z_optimal,inst):
source
MaximumWeightTwoStageSpanningTree.TwoStageSpanningTreeInstanceType
TwoStageSpanningTreeInstance

Contains all the relevant information defining a two stage spanning-tree instance

Fields

  • g::SimpleGraph{Int}: Graph
  • edge_index::SparseMatrixCSC{Int64, Int64}: edge_index[src(e),dst(e)] contains the Int index of edge e
  • nb_scenarios::Int:
  • first_stage_weights_matrix::SparseMatrixCSC{Float64, Int64}:
  • first_stage_weights_vector::Vector{Float64}:
  • second_stage_weights::Vector{SparseMatrixCSC{Float64, Int64}}:
source
MaximumWeightTwoStageSpanningTree.TwoStageSpanningTreeInstanceMethod
function TwoStageSpanningTreeInstance(g::AbstractGraph;nb_scenarios=1, first_max=10, second_max=20, negative_weights=false)

Constructor of TwoStageSpanningTreeInstance

  • g::AbstractGraph is the graph
  • first stage costs are randomly chosen between 1 and first_max
  • second stage costs are randomly chosen between 1 and second_max
  • weights are negative is negative_weights is true (multiplies weights above by -1)
source
MaximumWeightTwoStageSpanningTree.build_datasetMethod

function builddataset(; nbscenarios=5:5:20, firstmax=20:20, secondmax=10:5:30, seeds=1:5, gridsizes=4:6, solver=benderssolver, )

Build a training/valisation/test dataset for two stage spanning tree pipelines with a minimum weight spanning tree CO layer
source
MaximumWeightTwoStageSpanningTree.build_flow_graph_for_constraint_pricing!Function
function build_flow_graph_for_constraint_pricing!(
    g::MetaGraph,
    flow_graph::MetaDiGraph,
    base_graph_vertices_vertex_index_in_flow_graph = Dict(),
    base_graph_edges_vertex_index_in_flow_graph = Dict()
)

Constraint separated: ∑{e in E(X)}xe <= |X| - 1 for any ∅ ⊊ X ⊊ V`

  • g is the undirected graph on which the forest polytope is manipulated.
  • flow_graph is a MetaDiGraph. It should be empty in input. It is modified by the function and contains afterwards the MetaDiGraph used to separate the forest polytope constraint on f
source
MaximumWeightTwoStageSpanningTree.build_load_or_solveMethod
function build_load_or_solve(;
    graph=grid(5,5),
    seed=0,
    nb_scenarios=10,
    first_max=10,
    second_max=20,
    solver=lagrangian_heuristic_solver,
    load_and_save=true
)

Three solvers available

- `cut_solver`
- `benders_solver`
- `lagrangian_heuristic_solver`

return inst, val, sol with

- `inst`: the instance generated
- `val`: value of the solution computed
- `solution_computed`: solution computed
source
MaximumWeightTwoStageSpanningTree.lagrangian_dualMethod
lagrangian_dual(θ::AbstractArray; inst::TwoStageSpanningTreeInstance)

computes the value of the lagrangian dual function for TwoStageSpanningTreeInstance instance inst with duals θ

θ[edged_index(src(e),dst(e)),s] contains the value of the Lagrangian dual corresponding to edge e for scenario s

ChainRulesCore.rrule automatic differentiation with respect to θ works.

return (value of the solution computed),(edges in the solution computed)
source
MaximumWeightTwoStageSpanningTree.lagrangian_dual_stochasticMethod
lagrangian_dual(θ::AbstractArray; inst::TwoStageSpanningTreeInstance)

computes a stochastique (for one random scenario) value of the lagrangian dual function for TwoStageSpanningTreeInstance instance inst with duals θ

θ[edged_index(src(e),dst(e)),s] contains the value of the Lagrangian dual corresponding to edge e for scenario s

ChainRulesCore.rrule automatic differentiation with respect to θ works (and returns a stochastix gradient of the lagrangian dual function)

return (value of the solution computed),(edges in the solution computed)
source
MaximumWeightTwoStageSpanningTree.lagrangian_heuristicMethod
function lagrangian_heuristic(θ::AbstractArray; inst::TwoStageSpanningTreeInstance)

performs a lagrangian heuristic on TwoStageSpanningTree instance inst with duals θ.

θ[edge_index(src(e),dst(e)),s] contains the value of the Lagrangian dual corresponding to edge e for scenario s

return (value of the solution computed),(edges in the solution computed)
source
MaximumWeightTwoStageSpanningTree.maximum_weight_spanning_tree_single_scenario_layer_linear_encoderMethod
function maximum_weight_spanning_tree_single_scenario_layer_linear_encoder(inst::TwoStageSpanningTreeInstance)
    (inst::TwoStageSpanningTreeInstance)

returns X::Array{Float64} with `X[f,edge_index(inst,e),s]` containing the value of feature number `f` for edge `e` and scenario `s`

Features used: (all are homogeneous to a cost)
- first_stage_cost 
- second_stage_cost_quantile
- neighbors_first_stage_cost_quantile 
- neighbors_scenario_second_stage_cost_quantile 
- is_in_first_stage_x_first_stage_cost 
- is_in_second_stage_x_second_stage_cost_quantile 
- is_first_in_best_stage_x_best_stage_cost_quantile 
- is_second_in_best_stage_x_best_stage_cost_quantile 

For features with quantiles, the following quantiles are used: 0:0.1:1
source
MaximumWeightTwoStageSpanningTree.separate_forest_polytope_constraint_vertex_set_using_min_cut_MILP_formulation!Method
separate_forest_polytope_constraint_vertex_set_using_min_cut_MILP_formulation!(g::MetaGraph; MILP_solver=GLPK.Optimizer)

Uses a simple MILP to separate the constraints Constraint separated: ∑{e in E(X)} xe <= |X| - 1 for any ∅ ⊊ X ⊊ V`

  • g is a MetaGraph, :cb_val property contains the value of x_e

returns: found, X

  • found : boolean indicating if a violated constraint has been found
  • X : vertex set corresponding to the violated constraint
source
MaximumWeightTwoStageSpanningTree.separate_mst_Benders_cut!Method
separate_mst_Benders_cut!(g::MetaGraph ; MILP_solver=GLPK.Optimizer)

separates optimality and feasibility cuts for the Benders decomposition formulation of the MST problem

input: property :x_val contains the value of the master property :weight contains the value of the second stage

returns a boolean equals to true if the cut is a feasibility cut and false if it is an optimality cut (in that case, the cut might be satisfied)

The value of the duals μ and ν are stored in properties :mu and :nu of the g

source
MaximumWeightTwoStageSpanningTree.train_save_load_BBL_modelMethod

function trainsaveloadBBLmodel( ; traindataname, # Data name (to save on disk) traindata, # Data obtained with `buildorloadspanningtreeCOlayerdatasets()` nbDIRECTiterations=1000, # Nb iterations of the black box algorithm (DIRECT) perturbed=false, # Activate Gaussian perturbation nbperturbations=20, # Nb pertubation scenarios perturbationintensity=0.1, # Strength of the pertubation force_recompute=false # learn even if model available on disk )

source

Index