OptimalBranchingMIS

The maximum independent set problem

The maximum independent set problem is a classical combinatorial optimization problem, which is to find the largest subset of vertices in a graph such that no two vertices in the subset are adjacent.

Designing a branch-and-reduce algorithm using the optimal branching algorithm

To solve the MIS problem, we use the framework provided by the OptimalBranchingCore package to design a branch-and-reduce algorithm.

API

OptimalBranchingCore.apply_branchMethod
apply_branch(p::MISProblem{INTP, <:UnitWeight}, clause::Clause{INT}, vertices::Vector{T}) where {INTP<:Integer, INT<:Integer, T<:Integer}

Applies a clause to the given MISProblem and returns a new MISProblem with the left graph.

Arguments

  • p::MISProblem{INTP, <:UnitWeight}: The original problem.
  • clause::Clause{INT}: The clause to be applied.
  • vertices::Vector{T}: The vertices included in the clause.

Returns

  • A new MISProblem with the left graph.
source
OptimalBranchingCore.apply_branchMethod
apply_branch(p::MISProblem{INTP, <:Vector}, clause::Clause{INT}, vertices::Vector{T}) where {INTP<:Integer, INT<:Integer, T<:Integer}

Applies a clause to the given weighted MISProblem and returns a new weighted MISProblem with the left graph.

Arguments

  • p::MISProblem{INTP, <:Vector}: The original problem.
  • clause::Clause{INT}: The clause to be applied.
  • vertices::Vector{T}: The vertices included in the clause.

Returns

  • A new weighted MISProblem with the left graph.
source
OptimalBranchingCore.measureMethod
measure(p::MISProblem, ::D3Measure)

Calculates the D3 measure for the given MISProblem, which is defined as the sum of the maximum degree of each vertex minus 2, for all vertices in the graph.

Arguments

  • p::MISProblem: The problem instance containing the graph.

Returns

  • Int: The computed D3 measure value.
source
OptimalBranchingCore.measureMethod
measure(p::MISProblem, ::NumOfVertices)

Calculates the number of vertices in the given MISProblem.

Arguments

  • p::MISProblem: The problem instance containing the graph.

Returns

  • Int: The number of vertices in the graph.
source
OptimalBranchingCore.reduce_problemMethod
reduce_problem(::Type{R}, p::MISProblem, ::MISReducer) where R

Reduces the given MISProblem by removing vertices based on their degrees and returns a new MISProblem instance along with the count of removed vertices.

Arguments

  • p::MISProblem: The problem instance containing the graph to be reduced.
  • ::MISReducer: An implementation of the MISReducer abstract type.
  • ::Type{R}: The type of the result expected.

Returns

  • A tuple containing:
    • A new MISProblem instance with specified vertices removed.
    • An integer representing the count of removed vertices.

Description

The function checks the number of vertices in the graph:

  • If there are no vertices, it returns an empty instance and a count of 0.
  • If there is one vertex, it returns an empty instance and a count of 1 or weight of the vertex.
  • If there are two vertices, it returns an empty instance and a count based on the presence of an edge between them.
  • For graphs with more than two vertices, it calculates the degrees of the vertices and identifies the vertex with the minimum degree to determine which vertices to remove.
source
OptimalBranchingCore.size_reductionMethod
size_reduction(p::MISProblem{INT}, ::D3Measure, cl::Clause, variables::Vector) where {INT}

Calculates the size reduction for the given MISProblem using the D3 measure after applying the clause.

Arguments

  • p::MISProblem{INT}: The problem instance containing the graph.
  • cl::Clause: The clause to be applied.
  • variables::Vector: The variables included in the clause.

Returns

  • sum: The size reduction value.
source
OptimalBranchingCore.size_reductionMethod
size_reduction(p::MISProblem{INT}, ::NumOfVertices, cl::Clause, variables::Vector) where {INT}

Calculates the size reduction for the given MISProblem using the number of vertices after applying the clause.

Arguments

  • p::MISProblem{INT}: The problem instance containing the graph.
  • cl::Clause: The clause to be applied.
  • variables::Vector: The variables included in the clause.

Returns

  • sum: The size reduction value.
source
OptimalBranchingMIS.alphaMethod
alpha(g::SimpleGraph, weights::AbstractVector{WT}, openvertices::Vector{Int}) where WT

Compute the alpha tensor for a given weighted sub-graph.

Arguments

  • g::SimpleGraph: The input sub-graph.
  • weights::AbstractVector{WT}: The weights of the sub-graph.
  • openvertices::Vector{Int}: The open vertices of the sub-graph.

Returns

  • The alpha tensor.
source
OptimalBranchingMIS.clause_sizeMethod
clause_size(weights::Vector{WT}, bit_config, vertices::Vector) where WT

Compute the MIS size difference brought by the application of a clause for a given weighted graph.

Arguments

  • weights::Vector{WT}: The weights of the graph.
  • bit_config::Int: The bit configuration of the clause.
  • vertices::Vector{Int}: The vertices included in the clause.

Returns

  • The MIS size difference brought by the application of the clause.
source
OptimalBranchingMIS.closed_neighborsMethod
closed_neighbors(g::SimpleGraph, vertices::Vector{Int})

Returns a set of vertices that includes the input vertices as well as their open neighbors.

Arguments

  • g::SimpleGraph: The input graph.
  • vertices::Vector{Int}: The vertices for which closed neighbors are to be computed.

Returns

A set of vertices that includes the input vertices as well as their open neighbors.

source
OptimalBranchingMIS.counting_mis1Method
counting_mis1(g::SimpleGraph)

Calculates the maximum independent set (MIS) for a given SimpleGraph by first converting it into an EliminateGraph and then applying the counting_mis1 function.

Arguments

  • g::SimpleGraph: The input graph for which the maximum independent set is to be calculated.

Returns

  • MISCount: The size of the maximum independent set for the provided graph and the count of branches.
source
OptimalBranchingMIS.counting_xiao2013Method
counting_xiao2013(g::SimpleGraph)

This function counts the maximum independent set (MIS) in a given simple graph using the Xiao 2013 algorithm.

Arguments

  • g::SimpleGraph: A simple graph for which the maximum independent set is to be counted.

Returns

  • CountingMIS: An object representing the size of the maximum independent set and the count of branches.
source
OptimalBranchingMIS.counting_xiao2021Method
counting_xiao2021(g::SimpleGraph, weights::Vector{WT})

This function counts the maximum weighted independent set (MWIS) in a given simple graph with weights on vertices using the Xiao 2021 algorithm.

Arguments

  • g::SimpleGraph: A simple graph for which the maximum weighted independent set is to be counted.
  • weights::Vector{WT}: The weights of the vertices.

Returns

  • CountingMIS: An object representing the size of the maximum weighted independent set and the count of branches.
source
OptimalBranchingMIS.edge2vertexMethod
edge2vertex(g::SimpleGraph)

Connectivity between edges and vertices.

Arguments

  • g::SimpleGraph: The input graph.

Returns

  • A sparse matrix where the i-th row and j-th column is 1.0 if there is edge j had an end at vertex i.
source
OptimalBranchingMIS.graph_from_tuplesMethod
graph_from_tuples(n::Int, edgs)

Create a graph from a list of tuples.

Arguments

  • n::Int: The number of vertices.
  • edgs: The list of tuples.

Returns

  • The generated SimpleGraph g.
source
OptimalBranchingMIS.mis_branch_countMethod
mis_branch_count(g::AbstractGraph; weights::AbstractVector = UnitWeight(nv(g)), branching_strategy::BranchingStrategy = BranchingStrategy(table_solver = TensorNetworkSolver(), selector = MinBoundaryHighDegreeSelector(2, 6, 0), measure=D3Measure()), reducer=BasicReducer(), show_progress::Bool = false)

Calculate the size and the number of branches of the Maximum (Weighted) Independent Set (M(W)IS) for a given (vertex-weighted) graph.

Arguments

  • g::AbstractGraph: The graph for which the M(W)IS size and the number of branches are to be calculated.
  • weights::Vector: (optional) The weights of the vertices in the graph. It's set to be UnitWeight if the graph is not weighted. Defaults to UnitWeight(nv(g)).

Keyword Arguments

  • branching_strategy::BranchingStrategy: (optional) The branching strategy to be used. Defaults to a strategy using table_solver=TensorNetworkSolver, selector=MinBoundaryHighDegreeSelector(2, 6, 0), and measure=D3Measure.
  • reducer::AbstractReducer: (optional) The reducer to be applied. Defaults to BasicReducer.
  • show_progress::Bool: (optional) Whether to show the progress of the branching and reduction process. Defaults to false.

Returns

  • A tuple (size, count) where size is the size of the Maximum (Weighted) Independent Set and count is the number of branches.
source
OptimalBranchingMIS.mis_sizeMethod
mis_size(g::AbstractGraph; weights::AbstractVector = UnitWeight(nv(g)), branching_strategy::BranchingStrategy = BranchingStrategy(table_solver = TensorNetworkSolver(), selector = MinBoundaryHighDegreeSelector(2, 6, 0), measure=D3Measure()), reducer::AbstractReducer = BasicReducer(), show_progress::Bool = false)

Calculate the size of the Maximum (Weighted) Independent Set (M(W)IS) for a given (vertex-weighted) graph.

Arguments

  • g::AbstractGraph: The graph for which the M(W)IS size is to be calculated.
  • weights::Vector: (optional) The weights of the vertices in the graph. It's set to be UnitWeight if the graph is not weighted. Defaults to UnitWeight(nv(g)).

Keyword Arguments

  • branching_strategy::BranchingStrategy: (optional) The branching strategy to be used. Defaults to a strategy using table_solver=TensorNetworkSolver, selector=MinBoundaryHighDegreeSelector(2, 6, 0), and measure=D3Measure.
  • reducer::AbstractReducer: (optional) The reducer to be applied. Defaults to BasicReducer.
  • show_progress::Bool: (optional) Whether to show the progress of the branching and reduction process. Defaults to false.

Returns

  • An integer representing the size of the Maximum (Weighted) Independent Set for the given (vertex-weighted) graph.
source
OptimalBranchingMIS.neighbor_coverMethod
neighbor_cover(g::SimpleGraph, v::Int, k::Int)

Compute the neighbor cover of a vertex in a graph.

Arguments

  • g::SimpleGraph: The input graph.
  • v::Int: The vertex for which to compute the neighbor cover.
  • k::Int: The number of iterations to perform.

Returns

  • vertices: An array containing the vertices in the neighbor cover.
  • openvertices: An array containing the open vertices in the neighbor cover.
source
OptimalBranchingMIS.neighbors_2ndMethod
neighbors_2nd(g::SimpleGraph, v::Int)

Return the second-order neighbors of a vertex v in a simple graph g.

Arguments

  • g::SimpleGraph: The simple graph.
  • v::Int: The vertex.

Returns

  • Array{Int}: An array of second-order neighbors of v.
source
OptimalBranchingMIS.open_neighborsMethod
open_neighbors(g::SimpleGraph, vertices::Vector{Int})

Returns a vector of vertices in the graph g, which are neighbors of the given vertices and not in the given vertices.

Arguments

  • g::SimpleGraph: The graph in which to find the open neighbors.
  • vertices::Vector{Int}: The vertices for which to find the open neighbors.

Returns

A vector of open neighbors of the given vertices.

source
OptimalBranchingMIS.open_verticesMethod
open_vertices(g::SimpleGraph, vertices::Vector{Int})

Remove vertices from the given vector that are connected to all other vertices in the graph.

Arguments

  • g::SimpleGraph: The graph object.
  • vertices::Vector{Int}: The vector of vertices.

Returns

  • Vector{Int}: The open vertices.
source
OptimalBranchingMIS.reduce_graphMethod
reduce_graph(g::SimpleGraph, weights::UnitWeight, ::BasicReducer)

Reduce the graph into a simplified one by basic reduction rules used in classical algorithms like mis2.

Arguments

  • g::SimpleGraph: The input graph.
  • weights::UnitWeight: The weights of the vertices.
  • ::BasicReducer: The reducer to be applied.

Returns

  • A ReductionResult instance containing the reduced graph, weights, mis difference, and vertex map.
source
OptimalBranchingMIS.reduce_graphMethod
reduce_graph(g::SimpleGraph, weights::UnitWeight, ::XiaoReducer)

Reduce the graph into a simplified one by Xiao's reduction rules.

Arguments

  • g::SimpleGraph: The input graph.
  • weights::UnitWeight: The weights of the vertices.
  • ::XiaoReducer: The reducer to be applied.

Returns

  • A ReductionResult instance containing the reduced graph, weights, mis difference, and vertex map.
source
OptimalBranchingMIS.reduce_graphMethod
reduce_graph(g::SimpleGraph, weights::AbstractVector, ::NoReducer)

Reduce the graph into a simplified one by reduction rules. Actually does nothing, only for comparison with other reducers.

Arguments

  • g::SimpleGraph: The input graph.
  • weights::AbstractVector: The weights of the vertices.
  • ::NoReducer: The reducer to be applied.

Returns

  • A ReductionResult instance containing the reduced graph, weights, mis difference, and vertex map.
source
OptimalBranchingMIS.reduce_graphMethod
reduce_graph(g::SimpleGraph, weights::Vector{WT}, ::BasicReducer)

Reduce the graph into a simplified one by basic weighted_version reduction rules.

Arguments

  • g::SimpleGraph: The input graph.
  • weights::Vector{WT}: The weights of the vertices.
  • ::BasicReducer: The reducer to be applied.

Returns

  • A ReductionResult instance containing the reduced graph, weights, mis difference, and vertex map.
source
OptimalBranchingMIS.reduce_graphMethod
reduce_graph(g::SimpleGraph, weights::AbstractVector{WT}, ::XiaoReducer)

Reduce the graph into a simplified one by Xiao's reduction rules.

Arguments

  • g::SimpleGraph: The input graph.
  • weights::AbstractVector{WT}: The weights of the vertices.
  • ::XiaoReducer: The reducer to be applied.

Returns

  • A ReductionResult instance containing the reduced graph, weights, mis difference, and vertex map.
source
OptimalBranchingMIS.reduced_alphaMethod
reduced_alpha(g::SimpleGraph, weights::AbstractVector{WT}, openvertices::Vector{Int}) where WT

Compute the reduced alpha tensor for a given weighted sub-graph.

Arguments

  • g::SimpleGraph: The input sub-graph.
  • weights::AbstractVector{WT}: The weights of the sub-graph.
  • openvertices::Vector{Int}: The open vertices of the sub-graph.

Returns

  • The reduced alpha tensor.
source
OptimalBranchingMIS.reduced_alpha_configsMethod
reduced_alpha_configs(solver::TensorNetworkSolver, graph::SimpleGraph, weights::AbstractVector{WT}, openvertices::Vector{Int}, potential=nothing) where WT

Compute the truth table according to the non-zero entries of the reduced alpha tensor for a given weighted sub-graph.

Arguments

  • solver::TensorNetworkSolver: The solver to use.
  • graph::SimpleGraph: The input sub-graph.
  • weights::AbstractVector{WT}: The weights of the sub-graph.
  • openvertices::Vector{Int}: The open vertices of the sub-graph.
  • potential::Union{Nothing, Vector{WT}}: The potential of the open vertices, defined as the sum of the weights of the nearestneighbors of the open vertices.

Returns

  • The truth table.
source
OptimalBranchingMIS.removed_verticesMethod
removed_vertices(vertices::Vector{Int}, g::SimpleGraph, clause::Clause{N}) where N

Given a list of vertices, a graph, and a clause, this function returns a list of removed vertices.

The vertices argument is a vector of integers representing the vertices to consider. The g argument is a SimpleGraph object representing the graph. The clause argument is a Clause object representing a clause.

The function iterates over the vertices and checks if the corresponding bit in the clause.mask is 1. If it is, the vertex is added to the list of removed vertices (rvs). If the corresponding bit in the clause.val is also 1, the neighbors of the vertex are also added to rvs.

The function returns the list of removed vertices with duplicates removed.

source
OptimalBranchingMIS.D3MeasureType
D3Measure

A struct representing a measure that calculates the sum of the maximum degree minus 2 for each vertex in the graph.

Fields

  • None
source
OptimalBranchingMIS.MISCountType
struct MISCount

Represents the count of Maximum Independent Sets (MIS).

Fields

  • mis_size::Int: The size of the Maximum Independent Set.
  • mis_count::Int: The number of Maximum Independent Sets of that size.

Constructors

  • MISCount(mis_size::Int): Creates a MISCount with the given size and initializes the count to 1.
  • MISCount(mis_size::Int, mis_count::Int): Creates a MISCount with the specified size and count.
source
OptimalBranchingMIS.MISProblemType
mutable struct MISProblem{INT <: Integer, VT<:AbstractVector} <: AbstractProblem

Represents a Maximum (Weighted) Independent Set (M(W)IS) problem.

Fields

  • g::SimpleGraph: The graph associated with the M(W)IS problem.
  • weights::VT: The weights of the vertices in the graph. It's set to be UnitWeight if the graph is not weighted.

Methods

  • copy(p::MISProblem): Creates a copy of the given MISProblem.
  • Base.show(io::IO, p::MISProblem): Displays the number of vertices in the MISProblem.
source
OptimalBranchingMIS.MinBoundaryHighDegreeSelectorType
struct MinBoundaryHighDegreeSelector <: AbstractVertexSelector

The MinBoundaryHighDegreeSelector struct represents a strategy: - if exists a vertex with degree geq highdegreethreshold, then select it and its k-degree neighbors. - otherwise, select a subgraph with the minimum number of open vertices by k-layers of neighbors.

Fields

  • kb::Int: The number of layers of neighbors to consider when selecting the subgraph.
  • hd::Int: The threshold of degree for a vertex to be selected.
  • kd::Int: The number of layers of neighbors to consider when selecting the subgraph.
source
OptimalBranchingMIS.MinBoundarySelectorType
struct MinBoundarySelector <: AbstractSelector

The MinBoundarySelector struct represents a strategy for selecting a subgraph with the minimum number of open vertices by k-layers of neighbors.

Fields

  • k::Int: The number of layers of neighbors to consider when selecting the subgraph.
source
OptimalBranchingMIS.ReductionResultType
struct ReductionResult{VT<:AbstractVector, WT}

A struct to store the result of the reduction.

Fields

  • g::SimpleGraph{Int}: The reduced graph.
  • weights::VT: The weights of the reduced graph.
  • r::WT: The reduction value (MIS difference).
  • vmap::Vector{Int}: The vertex map from the vertices in the reduced graph to the vertices in the original graph.
source
OptimalBranchingMIS.SubsolverReducerType
struct SubsolverReducer{TR, TS} <: MISReducer

After the size of the problem is smaller than the threshold, use a subsolver to find reduction rules.

Fields

  • reducer::TR = XiaoReducer(): The reducer used to reduce the graph.
  • subsolver::Symbol = :xiao: subsolvers include :mis2, :xiao and :ip.
  • threshold::Int = 100: The threshold for using the subsolver.
source
OptimalBranchingMIS.TensorNetworkReducerType
struct TensorNetworkReducer <: MISReducer

A reducer that uses tensor network contraction to find reduction rules.

Fields

  • n_max::Int = 15: Maximum number of vertices to be included in the selected region
  • selector::Symbol = :mincut: Strategy for selecting vertices to contract. Options are:
    • :neighbor: Select vertices based on neighborhood
    • :mincut: Select vertices based on minimum cut provided by KaHyPar (as extension, requires using KaHyPar) !Note: KaHyPar may leads to memory leak!
  • measure::AbstractMeasure = NumOfVertices(): Measure used for kernelization
  • intersect_strategy::Symbol = :bfs: Strategy for intersecting clauses. Options are:
    • :bfs: Breadth-first search (gives the optimal result)
    • :dfs: Depth-first search (gives the first found non-zero intersection)
  • sub_reducer::MISReducer = BasicReducer(): Reducer applied to selected vertices before tensor network contraction
  • region_list::Dict{Int, Tuple{Vector{Int}, Vector{Int}}}: Store the selected region
  • recheck::Bool = false: Whether to recheck the vertices that have not been modified
  • folding::Bool = false: Whether to fold the vertices
source
OptimalBranchingMIS.TensorNetworkSolverType
TensorNetworkSolver
TensorNetworkSolver(; prune_by_env::Bool = true)

A struct representing a solver for tensor network problems. This struct serves as a specific implementation of the AbstractTableSolver type.

source