Reference

ProblemReductions.CircuitType
struct Circuit

A circuit expression is a sequence of assignments.

Fields

  • exprs::Vector{Assignment}: The assignments in the circuit.
source
ProblemReductions.CircuitSATType
struct CircuitSAT <: AbstractProblem

Circuit satisfiability problem, where the goal is to find an assignment that satisfies the circuit.

Fields

  • circuit::Circuit: The circuit expression in SSA form.
  • symbols::Vector{Symbol}: The variables in the circuit.
source
ProblemReductions.ColoringType
struct Coloring{K, WT<:(AbstractVector)} <: AbstractProblem
Coloring{K}(graph; weights=UnitWeight(nv(graph)))

The Vertex Coloring problem.

Positional arguments

  • graph is the problem graph.
  • weights are associated with the edges of the graph, default to UnitWeight(ne(graph)).
source
ProblemReductions.DominatingSetType
struct DominatingSet{GT<:Graphs.AbstractGraph} <: AbstractProblem
DominatingSet(graph; weights=UnitWeight())

The dominating set problem.

Positional arguments

  • graph is the problem graph.

We don't have weights for this problem.

source
ProblemReductions.FactoringType
Factoring <: AbstractProblem

Factoring problem. Given two numbers a and b, find the product a * b.

Fields

  • m::Int64

  • n::Int64

  • input::Int64

where m is the number of bits for the first number, n is the number of bits for the second number, and input is the number to factorize.

source
ProblemReductions.GridGraphType
struct GridGraph <: Graphs.AbstractGraph{Int64}

A grid graph is a graph in which the vertices are arranged in a grid and two vertices are connected by an edge if and only if they are adjacent in the grid.

Fields

  • grid::BitMatrix: a matrix of booleans, where true indicates the presence of an edge.
  • radius::Float64: the radius of the unit disk
source
ProblemReductions.HyperGraphType
struct HyperGraph <: Graphs.AbstractGraph{Int64}

A hypergraph is a generalization of a graph in which an edge can connect any number of vertices.

Fields

  • n::Int: the number of vertices
  • edges::Vector{Vector{Int}}: a vector of vectors of integers, where each vector represents a hyperedge connecting the vertices with the corresponding indices.
source
ProblemReductions.KSatisfiabilityType
struct KSatisfiability{K, T} <: ProblemReductions.AbstractSatisfiabilityProblem{T}

The satisfiability problem for k-SAT, where the goal is to find an assignment that satisfies the CNF.

Fields

  • variables::Vector{T}: The variables in the CNF.
  • cnf::CNF{T}: The CNF expression.
source
ProblemReductions.LogicGadgetType
struct LogicGadget{PT<:AbstractProblem}

The logic gadget defined on an computational model.

Fields

  • problem::PT: the computational model, e.g., SpinGlass.
  • inputs::Vector{Int}: the input variables.
  • outputs::Vector{Int}: the output variables.

References

source
ProblemReductions.MaxCutType
struct MaxCut{WT<:(AbstractVector)} <: AbstractProblem

The cutting problem.

  • In this problem, we would like to find the cut of the graph that maximizes the sum of the

weights of the edges that are cut.

Positional arguments

  • graph is the problem graph.
  • weights are associated with the edges of the graph. We have ensure that the weights are in the same order as the edges in edges(graph).
source
ProblemReductions.MaximalISType
struct MaximalIS{WT<:Union{UnitWeight, Vector}} <: AbstractProblem

The [maximal independent set]problem. In the constructor, weights are the weights of vertices.

Positional arguments

  • graph is the problem graph.
  • weights are associated with the vertices of the graph.
source
ProblemReductions.PaintShopType
struct PaintShop{LT} <: AbstractProblem

The binary paint shop problem.

Positional arguments

  • sequence is a vector of symbols, each symbol is associated with a color.
  • isfirst is a vector of boolean numbers, indicating whether the symbol is the first appearance in the sequence.
source
ProblemReductions.QUBOType
struct QUBO{T<:Real} <: AbstractProblem

The quadratic unconstrained binary optimization.

\[E = \sum_{i,j} Q_{ij} x_i x_j\]

where x_i \in \{0, 1\}.

Arguments

  • matrix::AbstractMatrix: the matrix Q of the QUBO problem.
source
ProblemReductions.ReductionCircuitToSpinGlassType
struct ReductionCircuitToSpinGlass{GT, T} <: AbstractReductionResult

The reduction result of a circuit to a spin glass problem.

Fields

  • num_source_vars::Int: the number of variables in the source circuit.
  • spinglass::SpinGlass{GT, T}: the spin glass problem.
  • variables::Vector{Int}: the variables in the spin glass problem.
source
ProblemReductions.ReductionFactoringToSatType
struct ReductionFactoringToSat <: AbstractReductionResult

The reduction result of a factoring problem to a CircuitSAT problem.

Fields

  • circuit::CircuitSAT: the CircuitSAT problem.
  • p::Vector{Int}: the first number to multiply (store bit locations)
  • q::Vector{Int}: the second number to multiply.
  • m::Vector{Int}: the result of the multiplication.
source
ProblemReductions.ReductionGraphType
ReductionGraph

A directed graph representing the reduction paths between different problems. A node represents a problem type, and an edge represents a reduction rule from one problem type to another.

Fields

  • graph::Graphs.SimpleGraphs.SimpleDiGraph{Int64}

  • nodes::Vector{Any}

  • method_table::Dict{Pair{Int64, Int64}, Method}

source
ProblemReductions.ReductionMaxCutToSpinGlassType
struct ReductionMaxCutToSpinGlass{GT, T} <: AbstractReductionResult

The reduction result of a maxcut to a spin glass problem.

Fields

  • spinglass::SpinGlass{GT, T}: the spin glass problem.

We only consider a simple reduction from MaxCut to SpinGlass(the graph must be SimpleGraph).

source
ProblemReductions.ReductionQUBOToSpinGlassType
struct ReductionQUBOToSpinGlass{GT, T} <: AbstractReductionResult

The reduction result of a qubo to a spin glass problem.

Fields

  • spinglass::SpinGlass{GT, T}: the spin glass problem.

We only consider a simple reduction from QUBO to SpinGlass(the graph must be SimpleGraph).

source
ProblemReductions.ReductionSATToDominatingSetType
struct ReductionSATToDominatingSet{GT<:Graphs.AbstractGraph} <: AbstractReductionResult

The reduction result of a general SAT problem to an Dominating Set problem.

Fields

  • target::DominatingSet

  • num_literals::Int64

  • num_clauses::Int64

source
ProblemReductions.ReductionSATToIndependentSetType
struct ReductionSATToIndependentSet{T, GT<:Graphs.AbstractGraph, WT<:(AbstractVector)} <: AbstractReductionResult

The reduction result of a general SAT problem to an Independent Set problem.

Fields

  • target::IndependentSet

  • literals::Array{BoolVar{T}, 1} where T

  • source_variables::Vector

  • num_clauses::Int64

source
ProblemReductions.ReductionSatToColoringType
struct ReductionSatToColoring{K, T, WT<:(AbstractVector)} <: AbstractReductionResult

The reduction result of a Sat problem to a Coloring problem.

Fields

  • Coloring{K, WT<:AbstractVector}: the coloring problem, where K is the number of colors and WT is the weights type.
  • varlabel, used to filter extra variables

Note: The coloring problem is a 3 coloring problem, in which a auxiliary color is used Auxiliary color => 2.

source
ProblemReductions.ReductionSpinGlassToMaxCutType
struct ReductionSpinGlassToMaxCut{WT} <: AbstractReductionResult

The reduction result of a spin glass to a maxcut problem.

Fields

  • maxcut::MaxCut{WT}: the MaxCut problem.
  • ancilla::Int: the ancilla vertex.
source
ProblemReductions.ReductionVertexCoveringToSetCoveringType
struct ReductionVertexCoveringToSetCovering{ET, WT<:(AbstractVector)} <: AbstractReductionResult

The reduction result of a vertex covering to a set covering problem.

Fields

  • setcovering::SetCovering{ET,WT}: the set covering problem, where ET is the sets type and WT is the weights type.
  • edgelabel: map each edge to a number in order to identify the edge (otherwise the vector would be cluttering)
source
ProblemReductions.SatisfiabilityType
struct Satisfiability{T} <: ProblemReductions.AbstractSatisfiabilityProblem{T}

The satisfiability problem.

Fields

  • cnf is a conjunctive normal form (CNF) for specifying the satisfiability problems.
  • weights are associated with clauses.
source
ProblemReductions.SetPackingType
struct SetPacking{ET} <: AbstractProblem

The set packing problem, a generalization of independent set problem to hypergraphs.

Positional arguments

  • elements is a vector of elements in the universe.
  • sets is a vector of vectors, each set is associated with a weight specified in weights.

Currently this type problem doesn't support weights.

source
ProblemReductions.SpinGlassType
struct SpinGlass{GT<:Graphs.AbstractGraph, WT<:(AbstractVector)} <: AbstractProblem
SpinGlass(graph::AbstractGraph, J, h=zeros(nv(graph)))

The spin-glass problem.

Positional arguments

  • graph is a graph object.
  • weights are associated with the edges.
source
ProblemReductions.StaticBitVectorType
StaticBitVector{N,C} = StaticElementVector{N,1,C}
StaticBitVector(x::AbstractVector)

Examples

julia> sb = StaticBitVector([1,0,0,1,1])
10011

julia> sb[3]
0x0000000000000000

julia> collect(Int, sb)
5-element Vector{Int64}:
 1
 0
 0
 1
 1
source
ProblemReductions.StaticElementVectorType
StaticElementVector{N,S,C}
StaticElementVector(nflavor::Int, x::AbstractVector)

N is the length of vector, C is the size of storage in unit of UInt64, S is the stride defined as the log2(# of flavors). When the number of flavors is 2, it is a StaticBitVector.

Fields

  • data is a tuple of UInt64 for storing the configuration of static elements.

Examples

julia> ev = StaticElementVector(3, [1,2,0,1,2])
12012

julia> ev[2]
0x0000000000000002

julia> collect(Int, ev)
5-element Vector{Int64}:
 1
 2
 0
 1
 2
source
ProblemReductions.TruthTableType
struct TruthTable{N, T}

The truth table.

Fields

  • inputs::Vector{T}: The input values.
  • outputs::Vector{T}: The output values.
  • values::Vector{BitStr{N, Int}}: The truth table values.

Examples

julia> tt = TruthTable(['a', 'b'], ['c'], [bit"0", bit"0", bit"0", bit"1"])
┌───┬───┬───┐
│ a │ b │ c │
├───┼───┼───┤
│ 0 │ 0 │ 0 │
│ 1 │ 0 │ 0 │
│ 0 │ 1 │ 0 │
│ 1 │ 1 │ 1 │
└───┴───┴───┘
source
ProblemReductions.UnitDiskGraphType
struct UnitDiskGraph{D, T} <: Graphs.AbstractGraph{Int64}

A unit disk graph is a graph in which the vertices are points in a plane and two vertices are connected by an edge if and only if the Euclidean distance between them is at most a given radius.

Fields

  • n::Int: the number of vertices
  • locations::Vector{NTuple{D, T}}: the locations of the vertices
  • radius::T: the radius of the unit disk
source
ProblemReductions.VertexCoveringType
struct VertexCovering{WT<:(AbstractVector)} <: AbstractProblem

Vertex covering is a problem that seeks to find a minimum set of vertices that cover all edges in a graph.

Positional arguments

  • graph is a graph object.
  • weights are associated with the vertices of the graph, default to UnitWeight(nv(graph)).
source
ProblemReductions.evaluateFunction
evaluate(problem::AbstractProblem, config) -> Real

Evaluate the energy of the problem given the configuration config. The lower the energy, the better the configuration.

source
ProblemReductions.evaluateMethod
evaluate(c::Coloring, config)

Compute the energy of the vertex coloring configuration config, the energy is the number of violated edges.

source
ProblemReductions.evaluateMethod
evaluate(c::IndependentSet, config)

Firstly, we count the edges connecting the input 'config' (a subset of vertices): If this number is zero, this 'config' corresponds to an Independent Set.

  • If the 'config' is an independent set, we return - (size(independent set));
  • If the 'config' is not an independent set, we return Inf.
source
ProblemReductions.evaluateMethod
evaluate(c::Matching, config)
Return Inf if the configuration is not a matching, otherwise return the sum of the weights of the edges in the matching.
source
ProblemReductions.evaluateMethod
evaluate(c::MaxCut, config)

Compute the cut weights for the vertex configuration config (an iterator). The energy is the sum of the weights of the edges that are cut.

source
ProblemReductions.evaluateMethod
evaluate(ps::PaintShop, config)

Returns the number of color switches. For example, if the sequence is abaccb ,there are three variables, then the config should be [1,0,1] or [0,1,0]. Here [1,0,1] means you want the first color for a and c is red, and the first color for b is blue.

source
ProblemReductions.evaluateMethod
evaluate(c::SetPacking, config)
  • First step: We check if config (a vector of boolean numbers as the mask of sets) is a set packing of sets;
  • Second step: If it is a set packing, we return - (size(set packing)); Otherwise, we return size(variables) + 1.
source
ProblemReductions.evaluateMethod
evaluate(c::VertexCovering, config)

return the weights of edge that is not covered but return typemax(eltype(weights)) if the edge is not covered. config is a vector of boolean numbers.

source
ProblemReductions.extract_multiple_solutionsMethod
extract_multiple_solutions(reduction::AbstractReductionResult, solution_set)

Extract multiple solutions together solution_set of the target problem to the original problem.

Arguments

  • reduction: The reduction result.
  • solution_set: The set of multiple solutions of the target problem.
source
ProblemReductions.extract_solutionFunction
extract_solution(reduction::AbstractReductionResult, solution)

Extract the solution solution of the target problem to the original problem.

Arguments

  • reduction: The reduction result.
  • solution: The solution of the target problem.
source
ProblemReductions.flavorsMethod
flavors(::Type{<:AbstractProblem}) -> Vector

Returns a vector of integers as the flavors (domain) of a degree of freedom.

source
ProblemReductions.implement_reduction_pathMethod
implement_reduction_path(rg::ReductionGraph, path::AbstractVector, problem::AbstractProblem)

Implement a reduction path on a problem. Returns a ConcatenatedReduction instance.

Arguments

  • path::AbstractVector: A sequence of problem types, each element is an AbstractProblem instance.
  • problem::AbstractProblem: The source problem, the type of which must be the same as the first node in the path.
source
ProblemReductions.is_matchingMethod
is_matching(graph::SimpleGraph, config)

Returns true if config is a valid matching on graph, and false if a vertex is double matched. config is a vector of boolean variables, which has one to one correspondence with edges(graph).

source
ProblemReductions.is_set_packingMethod
is_set_packing(sets::AbstractVector, config)

Return true if config (a vector of boolean numbers as the mask of sets) is a set packing of sets.

source
ProblemReductions.is_vertex_coloringMethod
is_vertex_coloring(graph::SimpleGraph, config)

Returns true if the coloring specified by config is a valid one, i.e. does not violate the contraints of vertices of an edges having different colors.

source
ProblemReductions.is_vertex_coveringMethod
is_vertex_covering(graph::SimpleGraph, config)

return true if the vertex configuration config is a vertex covering of the graph. Our judgement is based on the fact that for each edge, at least one of its vertices is selected.

source
ProblemReductions.onehotvMethod
onehotv(::Type{<:StaticElementVector}, i, v)
onehotv(::Type{<:StaticBitVector, i)

Returns a static element vector, with the value at location i being v (1 if not specified).

source
ProblemReductions.paint_shop_coloring_from_configMethod
paint_shop_coloring_from_config(p::PaintShop, config)

Returns a valid painting from the paint shop configuration (given by the configuration solvers). The config is a sequence of 0 and 1, where 0 means painting the first appearence of a car in red, and 1 means painting the first appearence of a car in blue.

source
ProblemReductions.reduce_sizeMethod
reduce_size(::Type{T}, ::Type{S}, size)

Return the size of the target problem T after reducing the source problem S to T.

Note

The problem size measure is problem dependent. Please check problem_size for the problem size measure.

Arguments

  • T: The target problem type.
  • S: The source problem type.
  • size: The size of the source problem.
source
ProblemReductions.reducetoFunction
reduceto(::Type{TA}, x::AbstractProblem)

Reduce the problem x to a target problem of type TA. Returns an instance of AbstractReductionResult.

Arguments

  • TA: The target problem type.
  • x: The original problem.
source
ProblemReductions.reduction_pathsMethod
reduction_paths([rg::ReductionGraph, ]S::Type, T::Type)

Find all reduction paths from problem type S to problem type T. Returns a list of paths, where each path is a sequence of problem types.

Arguments

  • rg::ReductionGraph: The reduction graph of type ReductionGraph.
  • S::Type: The source problem type.
  • T::Type: The target problem type.
source
ProblemReductions.set_parametersFunction
set_parameters(problem::AbstractProblem, parameters) -> AbstractProblem

Change the parameters for the problem and return a new problem instance.

source
ProblemReductions.spinglass_energyMethod
spinglass_energy(g::SimpleGraph, config; J, h)
spinglass_energy(cliques::AbstractVector{Vector{Int}}, config; weights)

Compute the spin glass state energy for the vertex configuration config. In the configuration, the spin ↑ is mapped to configuration 0, while spin ↓ is mapped to configuration 1. Let $G=(V,E)$ be the input graph, the hamiltonian is

\[H = \sum_{ij \in E} J_{ij} s_i s_j + \sum_{i \in V} h_i s_i,\]

where $s_i \in \{-1, 1\}$ stands for spin ↓ and spin ↑.

In the hypergraph case, the hamiltonian is

\[H = \sum_{c \in C} w_c \prod_{i \in c} s_i,\]

where $C$ is the set of cliques, and $w_c$ is the weight of the clique $c$.

source
ProblemReductions.spinglass_gadgetMethod
spinglass_gadget(::Val{:arraymul})

The array multiplier gadget.

    s_{i+1,j-1}  p_i
           \     |
        q_j ------------ q_j
                 |
    c_{i,j} ------------ c_{i-1,j}
                 |     \
                 p_i     s_{i,j} 
  • variables: pi, qj, pq, c{i-1,j}, s{i+1,j-1}, c{i,j}, s{i,j}
  • constraints: 2 * c{i,j} + s{i,j} = pi qj + c{i-1,j} + s{i+1,j-1}
source
ProblemReductions.truth_tableMethod
truth_table(ga::LogicGadget; variables=1:num_variables(ga.problem), solver=BruteForce())

Compute the truth table of a logic gadget.

Arguments

  • ga::LogicGadget: the logic gadget.

Keyword Arguments

  • variables::Vector{Int}: the variables to be displayed.
  • solver::AbstractSolver: the solver to be used.
source
ProblemReductions.variablesFunction
variables(problem::AbstractProblem) -> Vector

The degrees of freedoms in the computational problem. e.g. for the maximum independent set problems, they are the indices of vertices: 1, 2, 3..., while for the max cut problem, they are the edges.

source
ProblemReductions.@circuitMacro
@circuit circuit_expr

Construct a circuit expression from a block of assignments.

Examples

julia> @circuit begin
        x = a ∨ b
        y = x ∧ c
       end
Circuit:
| x = ∨(a, b)
| y = ∧(x, c)
source