Functions

KidneyExchange.Bellman_Ford_cycle_searchMethod
Bellman_Ford_cycle_search(
    graph,
    vertex_cost,
    source,
    K,
    is_covered,
    pred,
    d,
    constant_cost
)

Bellman-Ford style search for one positive cost cycle

Arguments

  • graph::SimpleDiGraph : The directed graph with cost on each arc
  • arc_cost::Matrix{Float64}:
  • source::Int : The local vertex index from which starts the search
  • K::Int : The maximal length of cycles
  • pred[k][v] contains u if u is a predecessor of v in a path of length k from source

Return values

  • cycle::Vector{Int}: the positive cycle found, [] if none
source
KidneyExchange.MIP_chain_searchFunction
MIP_chain_search(mip, graph, source, is_vertex, arc_cost)
MIP_chain_search(
    mip,
    graph,
    source,
    is_vertex,
    arc_cost,
    verbose
)

IP model with subtour elimination constraints. The constraints are the generalized cutset inequalities (GCS) and they are added dynamically in a row generation algorithm. Refer for instance to the following reference for a presentation of the GCS: (Taccari, 2016).

Arguments

  • mip::Model: The JuMP model initialized with the flow conservation constraints and the bound on the length of the chain
  • graph::SimpleDiGraph : The directed graph with cost on each arc
  • source::Int: Index of the source vertex
  • is_vertex::BitVector: For each vertex, indicates if it is in the considered subgraph
  • arc_cost::Vector{Float64}: Matrix of reduced costs of every arc

Return values

  • is_positivie_chain::Bool: True if a positive chain was found
  • chain::Vector{Int}: Positive chain that was found
source
KidneyExchange.add_column_to_masterMethod
add_column_to_master(column, mastermodel, tree_node)

Arguments

  • column::Column: Column to add to the master model
  • mastermodel::Model: JuMP model for current master problem
  • treenode::TreeNode: Information on current tree node

Return values: None

source
KidneyExchange.bfsFunction
bfs(g, source, K)
bfs(g, source, K, vertex_in_subgraph)
bfs(g, source, K, vertex_in_subgraph, d_from_vstar)

Breadth-first search to calculate the shortest path in terms of number of arcs from source to other vertices of graph g. The search stops when it is at level K even if there is vertex to be visited. Because we are trying to find cycles of at most K length, vertices of level >=K are never be considered in cycles.

Input parameters

  • g::AbstractGraph{T} : The graph
  • source::Int : The source vertex
  • K::Int: The maximum number of arcs
  • vertex_in_subgraph::Array{Bool}: Table indicating for each vertex if it is in the subgraph under consideration

Output parameters

  • dists::Vector{Float64}: shortest distance from source to vertices
source
KidneyExchange.bfs_reverseFunction
bfs_reverse(g, source, K)
bfs_reverse(g, source, K, vertex_in_subgraph)
bfs_reverse(g, source, K, vertex_in_subgraph, d_to_vstar)

Breadth-first search in the graph where arcs are reversed

Input parameters

  • g::AbstractGraph{T} : The graph
  • source::Int : The source vertex
  • K::Int: The maximum number of arcs
  • vertex_in_subgraph::Array{Bool}: Table indicating for each vertex if it is in the subgraph under consideration

Output parameters

  • dists::Vector{Float64}: shortest distance from source to vertices
source
KidneyExchange.branch_and_priceMethod
branch_and_price(
    instance,
    subgraphs,
    bp_params,
    timer,
    time_limit
)

Core function of the KEP solution with branch-and-price. It requires a parsed instance and the description of the graph copies. The column generation model is that of Riazcos-Alvarez et al (2020), but the many improvements have been added, in particular in the solution of the subproblem.

Arguments

  • instance::Instance: The parsed instance that is to be solved, it contains the KEP graph and the bounds on the length of covering cycles and chains
  • subgraphs::Graph_copies: Description of the graph copies of the extended edge formulation
  • bp_params::BP_params: solution parameters of the branch-and-price
  • timer::TimerOutput: a timer that will provide detail of where the computational time was spent during the branch-and-price
  • time_limit::Float64: time limit of the algorithm, including parsing and prepreprocessing times

#Output parametes

  • bp_status::BP_status: Structure containing every relevant information on the execution of the algorithm (including the optimal solution)
source
KidneyExchange.branch_on_arcFunction
branch_on_arc(
    arc_to_branch,
    master,
    is_cg_branching,
    tree,
    current_node,
    column_pool,
    node_count
)
branch_on_arc(
    arc_to_branch,
    master,
    is_cg_branching,
    tree,
    current_node,
    column_pool,
    node_count,
    verbose
)

Update the branch-and-bound tree with two new nodes by branching on the given arc with the specified branching (only on master problem or both in master and in subproblem)

Input parameters

  • arc_to_branch::Pair{Int,Int} : The arc to be branched
  • master::Model: Master model where the branching constraints are added
  • is_cg_branching::Bool: True if the branching impacts the subproblem of the coluln generation, false if it impacts only the master problem
  • tree::Vector{TreeNode}: Branch-and-bound tree
  • current_node::TreeNode: Branch-and-bound node currently treated
  • column_pool::Vector{Column}: Pool of all columns in current master problem
  • node_count::Int: Number of BP nodes enumerated until now
source
KidneyExchange.branch_on_vertexFunction
branch_on_vertex(
    vertex_to_branch,
    master,
    tree,
    current_node,
    column_pool,
    node_count
)
branch_on_vertex(
    vertex_to_branch,
    master,
    tree,
    current_node,
    column_pool,
    node_count,
    verbose
)

Update the branch-and-bound tree with two new nodes by branching on the given vertex.

Input parameters

  • vertex_to_branch::Int : The arc to be branched
  • master::Model: Master model where the branching constraints are added
  • tree::Vector{TreeNode}: Branch-and-bound tree
  • current_node::TreeNode: Branch-and-bound node currently treated
  • column_pool::Vector{Column}: Pool of all columns in current master problem

** node_count::Int: Number of BP nodes enumerated until now

source
KidneyExchange.build_hpief_mipFunction
build_hpief_mip(instance, subgraphs, params)
build_hpief_mip(instance, subgraphs, params, maxtime)

A compact MIP formulation originally proposed by (Dickerson et al., 2016). The hybrid position indexed edge formulation (HPIEF) handles cycles and chains in the same formulation.

source
KidneyExchange.build_reduced_extended_edge_mipFunction
build_reduced_extended_edge_mip(instance, subgraphs, params)
build_reduced_extended_edge_mip(
    instance,
    subgraphs,
    params,
    maxtime
)

A compact MIP formulation originally proposed by (Constantino et al., 2013) for the cycles-only variant of the problem. Here it is adapted to the graph copies based on FVS and chains are considered with position-indexed variables.

source
KidneyExchange.canGiveToMethod
canGiveTo(donorBT, patientBT)

Can a person with blood type donorBT give a kidney to a patient with blood type patientBT?

source
KidneyExchange.check_used_verticesMethod
check_used_vertices(cycle, cycles)

Gives a set of indices of cycle in cycles who have common vertices with cycle

Arguments

  • cycle::Array{Int,1}: The vertices of cycle
  • cycles::Array{Array{Int,1},1}: The set of cycles

#Output parametes

  • redundant_cycle_indices::Array{Int,1}: the set of indices of cycle in cycles who have common vertices with cycle
source
KidneyExchange.compute_arc_flowMethod
compute_arc_flow(mastersol, node_columns)

Calculates the relaxed value of the decision variable x(i,j) of arc (i->j) in A. x(i,j) whose value >0 is stored in x::Dict{Pair{Int,Int}, Float64}

Arguments

  • mastersol::Vector{Float64}: solution value of the master problem, ie: y[c] for c in cycles
  • node_columns::Vector{Column}: The corresponding cycles of current node

Return values

  • x:: Dict{Pair{Int,Int}, Float64}: value of x_(i,j) of arc (i->j)
source
KidneyExchange.create_chain_mipFunction
create_chain_mip(graph, L, optimizer, time_limit)
create_chain_mip(
    graph,
    L,
    optimizer,
    time_limit,
    nb_threads
)

Initialize the MIP model with cycle constraint generation for the optimal search of positive cost chains. Only one model is created for every copy to save a great amount of initialization time and memory. The model will then need to be modified for each graph copy to keep only the vertices of the graph and select the right source vertex

Arguments

  • graph::SimpleDiGraph : The directed graph with cost on each arc
  • L::Int: The maximal length of chains
  • optimizer::String: Name of the MIP sover
  • time_limit::Real: Time limit of the solver
  • gurobi_env: Gurobi environment if Gurobi is used (avoids many messages from the solver)

Return values

  • mip::Model: Initial JuMP model for the search of a positive chain
source
KidneyExchange.create_modelFunction
create_model(time_limit, optimizer)
create_model(time_limit, optimizer, is_integer)
create_model(time_limit, optimizer, is_integer, verbose)
create_model(
    time_limit,
    optimizer,
    is_integer,
    verbose,
    nb_threads
)
source
KidneyExchange.generateAltruistMethod
generateAltruist(pool_gen)

Random rolls an altruistic donor (donor with no attached patient)

  • ID unique identifier for the vertex
  • return altruistic donor vertex KPDVertexAltruist
source
KidneyExchange.generatePairMethod
generatePair(pool_gen)

Randomly rolls a patient-donor pair (possibly compatible or incompatible)

  • ID unique identifier for the vertex
  • return a patient-donor pair KPDVertexPair
source
KidneyExchange.generatePraIncompatibilityMethod
generatePraIncompatibility(pool_gen, isWifePatient)

Randomly generates CPRA (Calculated Panel Reactive Antibody) for a patient-donor pair, using the Saidman method. If the patient is the donor's wife, then CPRA is increased.

  • isWifePatient is the patent the wife of the donor?
  • return scaled CPRA double value between 0 and 1.0
source
KidneyExchange.generate_complete_benchmarkMethod
generate_complete_benchmark()

Generate all the instances that are used in addition to the PrefLib to assess the solution methods of this package in the article describing the package (Jérémy Omer, Ayse N Arslan, Fulin Yan. KidneyExchange.jl: A Julia package for solving the kidney exchange problem with branch-and-price. 2022. ⟨hal-03830810⟩).

The instances are stored in the subdirectories heterogeneous, sparse and saidman of the data directory. Beware that if using the package with Pkg.add(), the data directory is in the directory where the package is stored.

source
KidneyExchange.generate_heterogeneous_instanceMethod
generate_heterogeneous_instance(
    nb_pairs,
    nb_altruists,
    index
)

Generate a KEP dataset and write in two text files with extensions .dat and .wmd; the graph generator is based on the following article: Kidney Exchange in Dynamic Sparse Heterogeneous Pools: (Ashlagi et al., 2013).

Arguments

  • nb_pairs::Int: number of pairs of incompatible donor and patient
  • nb_altruists::Int: number of altruist donors
  • index::Int: index of the dataset, which will appear in the filename
source
KidneyExchange.generate_heterogeneous_kep_graphFunction
generate_heterogeneous_kep_graph(nb_pairs, nb_altruists)
generate_heterogeneous_kep_graph(
    nb_pairs,
    nb_altruists,
    pct_easy_to_match
)

Compatibility graph generator based on the following paper: Kidney Exchange in Dynamic Sparse Heterogeneous Pools. Itai Ashlagi, Patrick Jaillet, Vahideh H. Manshadi. EC-2013. (Extended abstract.) @author John P. Dickerson

source
KidneyExchange.generate_saidman_instanceMethod
generate_saidman_instance(nb_pairs, nb_altruists, index)

Generate a KEP dataset and write in two text files with extensions .dat and .wmd. This generator is usually refererred to as the "Saidman Generator".

Arguments

  • nb_pairs::Int: number of pairs of incompatible donor and patient
  • nb_altruists::Int: number of altruist donors
  • index::Int: index of the dataset, which will appear in the filename

Reference

(Saidman et al., Sep 2014)

source
KidneyExchange.generate_sparse_unos_instanceMethod
generate_sparse_unos_instance(nb_pairs, nb_altruists, index)

Generate a KEP dataset and write in two text files with extensions .dat and .wmd. This generator is obtained by applying the Saidman generator with small modifications in the probabilities of compatibility to roughtly mimic the UNOS pool as of April 15, 2013. The corresponding data was taken from the KPD Work Group Data Analysis - CMR - June 2013 report. Seel also (Dickerson et al., 2012) for the motivation and description of this generator.

Arguments

  • nb_pairs::Int: number of pairs of incompatible donor and patient
  • nb_altruists::Int: number of altruist donors
  • index::Int: index of the dataset, which will appear in the filename
source
KidneyExchange.get_branching_arcMethod
get_branching_arc(column_flow, pief_flow)

Find a fractional arc to branch. The fractional arc closest to 0.5 will be selected to branch, if there is no fractional arc in the solution

Input parameters

  • column_flow::Dict{Pair{Int,Int}, Float64} : The dictionary containing the nonzero flow on each arc due to the selection of columns
  • pief_flow::Dict{Pair{Int,Int}, Float64} : The dictionary containing the flow on each arc from the pief model for chain search (if applicable)

Output parameters

  • arc_to_branch::Pair{Int,Int} : The arc to be branched
  • is_cg_branching::Bool: True if the branching impacts the subproblem of the coluln generation, false if it impacts only the master problem
source
KidneyExchange.get_branching_vertexMethod
get_branching_vertex(graph, column_flow, pief_flow)

Find a fractional vertex cover to branch. The fractional vertex closest to 0.5 will be selected to branch

Input parameters

  • column_flow::Dict{Pair{Int,Int}, Float64} : The dictionary containing the nonzero flow on each arc due to the selection of columns
  • pief_flow::Dict{Pair{Int,Int}, Float64} : The dictionary containing the flow on each arc from the pief model for chain search (if applicable)

Output parameters

  • vertex_to_branch::Int : The vertex to be branched on, nothing if none was found
source
KidneyExchange.get_feasible_solutionMethod
get_feasible_solution(fractional_solution, node_columns)

Extract a integer feasible solution from the fractional solution by conserving a set of vertex-disjoint cycles or chains

Arguments

  • fractional_solution::Array{Array{Float64,1},1}: The set of value of variables indicating whether or not cycles and chains are selected
  • node_columns::Array{Array{Column,1},1}: The set of cycles and chains associated with variables

Output parameters

  • feas_val::Real: The objective value of the feasible solution
  • columns::Array{Array{Int,1},1}: The set of selected cycles and chains of the feasible solution
source
KidneyExchange.initialize_column_poolFunction
initialize_column_pool(instance, column_pool)
initialize_column_pool(
    instance,
    column_pool,
    max_cycle_length
)
initialize_column_pool(
    instance,
    column_pool,
    max_cycle_length,
    max_nb_cols
)

Initialize the pool of columns for the branch-and-price by enumerating all k-cycles up to an input maximum value of k.

Input parameters

*instance::Instance: The instance *column_pool::Vector{Column}: Pool of columns where initial columns will be pushed *max_cycle_length::Int: Maximum length of the cycles to be enumerated at initialization (0 if initialization is deactivated)

  • max_nb_cols::Int : Maximum number of columns added for each vertex, 0 if all columns are added

Output parameters

  • column_pool::Vector{Column} : set of columns initially added to the master problem
source
KidneyExchange.initialize_master_IPFunction
initialize_master_IP(instance, column_pool)
initialize_master_IP(instance, column_pool, bp_params)
initialize_master_IP(
    instance,
    column_pool,
    bp_params,
    time_limit
)

Initialization of the master problem with integer columns

Arguments

  • instance::Instance: The parsed instance that is to be solved, it contains the KEP graph and the bounds on the length of covering cycles and chains.
  • column_pool::Vector{Column}: the set of initial columns of the master
  • bp_params::BP_params: parameters of the branch-and-price
  • time_limit::Float64: time limit for each solution of the master relaxation

Return values

  • master::Model: the model of the restricted master problem
source
KidneyExchange.isDonorSpouseMethod
isDonorSpouse(pool_gen)
  • Draws a random spousal relationship between donor and patient
  • return true if willing donor is patient's spouse, false otherwise
source
KidneyExchange.isPositiveCrossmatchMethod
isPositiveCrossmatch(pr_PraIncompatibility)

Random roll to see if a patient and donor are crossmatch compatible

  • pr_PraIncompatibility: probability of a PRA-based incompatibility
  • return true is simulated positive crossmatch, false otherwise
source
KidneyExchange.node_masterFunction
node_master(instance, column_pool)
node_master(instance, column_pool, bp_params)
node_master(instance, column_pool, bp_params, time_limit)

Initialization of the restricted master problem

Arguments

  • instance::Instance: The parsed instance that is to be solved, it contains the KEP graph and the bounds on the length of covering cycles and chains.
  • column_pool::Vector{Column}: the set of initial columns of the master
  • bp_params::BP_params: parameters of the branch-and-price
  • time_limit::Float64: time limit for each solution of the master relaxation

Return values

  • master::Model: the model of the restricted master problem
source
KidneyExchange.preprocess_graph_copiesFunction
preprocess_graph_copies(instance)
preprocess_graph_copies(instance, reduce_arcs)
preprocess_graph_copies(
    instance,
    reduce_arcs,
    reduce_vertices
)
preprocess_graph_copies(
    instance,
    reduce_arcs,
    reduce_vertices,
    fvs
)

Get the data of each graph copy used in the subproblems of the column generation. The graph copies are charaterized by a source vertex and the set of vertices (optionally set of arcs) kept in the copy. Each altruist donor is the source of one graph copy. Preprocessing the vertices that can be reached with a given maximum number of arcs does not bring any improvement to the subproblem solution, so we do nothing at this stage.

For the other copies, either each donor/patient pair is the source of one copy, or we compute a feedback vertex set (FVS) that heuristically minimizes the cardinality of the set. Each vertex of the FVS will then be the source of one copy in the column generation subproblems. In both cases, we then compute distances from and to the source in each subgraph. These will be used in the Bellman-Ford search at each solution of the subproblems.

Input parameters

  • instance::Instance : The structure describing the characteristics of an instance including the compatibility graph
  • reduce_arcs::Bool : A boolean parameter indicating whether or not arcs that are not useful in a copy should be eliminated from it (default value is false)
  • reduce_vertices::Bool : A boolean parameter indicating whether or not vertices that are not useful in copy should be eliminated from it (default value is true)
  • fvs::Bool : True if the copies correspond to an FVS, false if there is one copy per donor/patient pair

Output parameters

  • subgraphs::Graph_copies : The structure describing the preprocessed data of the graph copies
source
KidneyExchange.print_and_check_solutionFunction
print_and_check_solution(cycles, chains, instance)
print_and_check_solution(cycles, chains, instance, verbose)

Function that prints the main characteristics of a solution and check that it satisfies all teh constraints

Input parameters

  • cycles::Vector{Vector{Int}}: list of cycles in the solution to display
  • chains::Vector{Vector{Int}}: list of chains in the solution to display
  • instance::Instance: KEP instance whose solution it is
  • verbose::Bool: true if the characteristics of the solution are printed, false if only the verification need be done

Output parameters

None

source
KidneyExchange.process_nodeMethod
process_node(
    tree_node,
    instance,
    mastermodel,
    subgraphs,
    bp_status,
    column_pool,
    bp_params,
    master_IP,
    timer,
    time_limit
)

Columns generation of the node

Input parameter

  • tree_node::TreeNode: Branch-and-bound node to solve with column generation

Output parameter

  • column_flow::Dict{Pair{Int,Int}, Float64}: value of x_(i,j) of arc (i->j)`
source
KidneyExchange.read_kep_fileMethod
read_kep_file(wmd_file, dat_file)

Contruct a KEP graph from a .wmd and a .dat input files. The format of these files can be found on the PrefLib website (see here for more details). The .dat file provides individual information on the patient and donor of each pair such as blood type. The .wmd file describes the edges of the KEP graph. For example, in the first instance of the benchmark, the .dat looks like this :

	Pair,Patient,Donor,Wife-P?,%Pra,Out-Deg,Altruist
	1,A,B,0,0.05,2,0
	2,O,A,0,0.05,4,0
	3,A,B,0,0.05,2,0
	4,O,A,1,0.5875,3,0
	5,B,AB,0,0.05,0,0
	6,B,A,0,0.45,3,0
	7,O,A,0,0.45,4,0
	8,B,A,0,0.45,3,0
	9,O,A,0,0.45,4,0
	10,O,O,1,0.2875,11,0
	11,A,AB,0,0.05,0,0
	12,A,A,1,0.5875,3,0
	13,O,O,1,0.925,10,0
	14,O,A,0,0.45,4,0
	15,O,A,0,0.05,3,0
	16,O,B,1,0.2875,3,0

And a preview of the .wmd file (including only a subset of the arcs) looks like this:

	# FILE NAME: 00036-00000001.wmd
	# TITLE: Kidney Matching - 16 with 0
	# DESCRIPTION:
	# DATA TYPE: wmd
	# MODIFICATION TYPE: synthetic
	# RELATES TO:
	# RELATED FILES: 00036-00000001.dat
	# PUBLICATION DATE: 2013-08-17
	# MODIFICATION DATE: 2022-09-16
	# NUMBER ALTERNATIVES: 16
	# NUMBER EDGES: 59
	# ALTERNATIVE NAME 1: Pair 1
	# ALTERNATIVE NAME 2: Pair 2
	# ALTERNATIVE NAME 3: Pair 3
	# ALTERNATIVE NAME 4: Pair 4
	...
	1, 5, 1.0
	1, 6, 1.0
	2, 1, 1.0
	2, 3, 1.0
	2, 11, 1.0
	2, 12, 1.0
	3, 5, 1.0
	3, 8, 1.0
	4, 1, 1.0
	4, 3, 1.0
	...

Parameters

  • wmd_file::String : Absolute path of the .wmd file.
  • dat_file::String : Absolute path of the .dat file.
source
KidneyExchange.read_wmd_fileMethod
read_wmd_file(filepath)

Contruct a KEP graph from a .wmd file. The format of these files can be found on the PrefLib website. The extra .dat file mentionned in .wmd file metadata must be provided. It contains individual information on the patient and donor of each pair such as blood type. The .wmd file describes the edges of the KEP graph. For example, in the first instance of the benchmark, the .dat looks like this :

	Pair,Patient,Donor,Wife-P?,%Pra,Out-Deg,Altruist
	1,A,B,0,0.05,2,0
	2,O,A,0,0.05,4,0
	3,A,B,0,0.05,2,0
	4,O,A,1,0.5875,3,0
	5,B,AB,0,0.05,0,0
	6,B,A,0,0.45,3,0
	7,O,A,0,0.45,4,0
	8,B,A,0,0.45,3,0
	9,O,A,0,0.45,4,0
	10,O,O,1,0.2875,11,0
	11,A,AB,0,0.05,0,0
	12,A,A,1,0.5875,3,0
	13,O,O,1,0.925,10,0
	14,O,A,0,0.45,4,0
	15,O,A,0,0.05,3,0
	16,O,B,1,0.2875,3,0

And a preview of the .wmd file (including only a subset of the arcs) looks like this:

	# FILE NAME: 00036-00000001.wmd
	# TITLE: Kidney Matching - 16 with 0
	# DESCRIPTION:
	# DATA TYPE: wmd
	# MODIFICATION TYPE: synthetic
	# RELATES TO:
	# RELATED FILES: 00036-00000001.dat
	# PUBLICATION DATE: 2013-08-17
	# MODIFICATION DATE: 2022-09-16
	# NUMBER ALTERNATIVES: 16
	# NUMBER EDGES: 59
	# ALTERNATIVE NAME 1: Pair 1
	# ALTERNATIVE NAME 2: Pair 2
	# ALTERNATIVE NAME 3: Pair 3
	# ALTERNATIVE NAME 4: Pair 4
	...
	1, 5, 1.0
	1, 6, 1.0
	2, 1, 1.0
	2, 3, 1.0
	2, 11, 1.0
	2, 12, 1.0
	3, 5, 1.0
	3, 8, 1.0
	4, 1, 1.0
	4, 3, 1.0
	...

Parameters

  • wmd_file::String : Absolute path of the .wmd file.
source
KidneyExchange.relaxed_arcFunction
relaxed_arc(instance, params)
relaxed_arc(instance, params, maxtime)

Deterministic relaxed-arc formulation for a KEP problem. No cycle length limitation is imposed.

Parameters

  • instance::Instance: KEP instance to solve
  • maxtime::Real=60 : Maximum solving time in seconds
source
KidneyExchange.solve_with_BPFunction
solve_with_BP(filename, K, L)
solve_with_BP(filename, K, L, bp_params)
solve_with_BP(filename, K, L, bp_params, timer)
solve_with_BP(filename, K, L, bp_params, timer, time_limit)

This is the main function to call for an execution of the branch-and-price algorithm on input file with given bounds on the length of covering cycles and chains and given options

Arguments

  • filename::String: path of the input data files, this should include the name of the files, but not the .dat and .wmd extensions
  • K::Int: maximum length of cycles
  • L::Int: maximum length of chains
  • bp_params::BP_params: solution parameters of the branch-and-price
  • timer::TimerOutput: a timer that will provide detail of where the computational time was spent during the branch-and-price
  • time_limit::Float64: time limit of the algorithm, including parsing and prepreprocessing times

#Output parametes

  • instance::Instance: The parsed instance that is to be solved, it contains the KEP graph and the bounds on the length of covering cycles and chains.
  • subgraphs::Graph_copies: Description of the graph copies of the extended edge formulation
  • bp_status::BP_status: Structure containing every relevant information on the execution of the algorithm (including the optimal solution)
source
KidneyExchange.solve_with_mipFunction
solve_with_mip(filename, K, L)
solve_with_mip(filename, K, L, params)
solve_with_mip(filename, K, L, params, timer)
solve_with_mip(filename, K, L, params, timer, time_limit)

This is the main function to call to solve the input instance with given bounds on the length of covering cycles and chains and given options

Arguments

  • filename::String: path of the input data files, this should include the name of the files, but not the .dat and .wmd extensions
  • K::Int: maximum length of cycles
  • L::Int: Maximum length of chains
  • params::MIP_params: parameters of the MIP model and solver
  • timer::TimerOutput: a timer that will provide detail of where the computational time was spent during the branch-and-price
  • time_limit::Float64: time limit of the algorithm, including parsing and prepreprocessing times

Output parametes

  • instance::Instance: The parsed instance that is to be solved, it contains the KEP graph and the bounds on the length of covering cycles and chains.
  • subgraphs::Graph_copies: Description of the graph copies of the extended edge formulation
  • bp_status::BP_status: Structure containing every relevant information on the execution of the algorithm (including the optimal solution)
source
KidneyExchange.traverse_predsMethod
traverse_preds(v, pred, n)

Inner function for BellmanFordcycle_search to return path of exactly (n-1) of length ending at v from the table of predecessor

Arguments

  • v::Int : The ending vertex
  • pred::Array{Array{Tuple{Int,Int},1},1}: table of predecessor containing tuples of predecessor and legnth of path to arrive the vertex
  • n::Int: the length of the path

Return values

  • c=Array{Int,1}: path of exactly (n-1) of length ending at v
source
KidneyExchange.write_kep_fileMethod
write_kep_file(
    kep_graph,
    edge_weight,
    donorBT,
    patientBT,
    wifeP,
    patientPRA,
    is_altruist,
    file_name
)

Write a .wmd and a .dat files to store the input KEP graph in the same form as in Preflib. Store the files in ./data directory.

Parameters

  • kep_graph::SimpleDiGraph: KEP graph to write
  • edge_weight::Matrix{Float64}: weight of each eadge of the KEP graph
  • donorBT::Vector{Blood_type}: blood type of the donor of each pair
  • patientBT::Vector{Blood_type}: blood type of the patient of each pair
  • wifeP::BitArray: for each pair true iff the donor and patient are married
  • patientPRA::Vector{Float64}: PRA of the patient of each pair
  • is_altruist::BitArray: for each pair, true if the donor is an altruist donor
  • file_name::String: Name of the files (before the extensions)
source
KidneyExchange.write_wmd_fileMethod
write_wmd_file(
    graph,
    weights,
    is_altruist,
    file_name,
    title,
    description
)

Write a .wmd and a .dat files to store the input KEP graph in the same form as in Preflib.

source