Reference

ProblemReductions.CircuitType
struct Circuit

A circuit expression is a sequence of assignments.

Fields

  • exprs::Vector{Assignment}: The assignments in the circuit.
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.LocalConstraintType
struct LocalConstraint

A constraint for specifying a ConstraintSatisfactionProblem, which is defined on finite domain variables.

Fields

  • num_flavors: the number of flavors (domain) of a degree of freedom.
  • variables: the indices of the variables involved in the constraint.
  • specification: a boolean vector of length num_flavors^length(variables), specifying whether a configuration is valid.
  • strides: the strides of the variables, to index the specification vector.
source
ProblemReductions.LocalSolutionSizeType
struct LocalSolutionSize{T}

Problem size defined on a subset of variables of a ConstraintSatisfactionProblem.

Fields

  • num_flavors: the number of flavors (domain) of a degree of freedom.
  • variables: the indices of the variables involved in the constraint.
  • specification: a vector of size num_flavors^length(variables), specifying the local solution sizes.
  • strides: the strides of the variables, to index the specification vector.
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.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
struct 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}, Function}

source
ProblemReductions.ReductionMatchingToSetPackingType
struct ReductionMatchingToSetPacking{ET, T, WT<:AbstractArray{T, 1}} <: AbstractReductionResult

The reduction result of a vertex matching to a set packing problem.

Fields

  • SetPacking{WT<:AbstractVector{Int}}: the target set packing problem
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.ReductionPathType
struct ReductionPath

A sequence of reductions.

Fields

  • nodes::Vector{AbstractProblem}: The sequence of problem types.
  • methods::Vector{Function}: The sequence of methods used to reduce the problems.
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.ReductionSATToCircuitType
struct ReductionSATToCircuit <: AbstractReductionResult

The reduction result of an SAT problem o a Circuit SAT problem.

Fields

  • target::CircuitSAT: the target problem.

    • target::CircuitSAT

    • sat_symbols::Vector{Symbol}

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{S, GT<:Graphs.AbstractGraph, T, WT<:AbstractArray{T, 1}} <: AbstractReductionResult

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

Fields

  • target::IndependentSet

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

  • source_variables::Vector

  • num_clauses::Int64

source
ProblemReductions.ReductionSatToColoringType
struct ReductionSatToColoring{K, T, WT<:AbstractArray{T, 1}} <: AbstractReductionResult

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

Fields

  • Coloring{K, T, WT<:AbstractVector{T}}: the coloring problem, where K is the number of colors and WT is the weights type.
  • posvertices, a map from literal to vertex index
  • negvertices, a map from negative literal to vertex index

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

source
ProblemReductions.ReductionSetPackingToIndependentSetType
struct ReductionSetPackingToIndependentSet{GT, ET} <: AbstractReductionResult

The reduction result of a Set Packing problem to an Independent Set problem.

Fields

  • target::IndependentSet{GT, T} where {GT, T}

  • subset_list::Array{Vector{ET}, 1} where ET

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, T, WT<:AbstractArray{T, 1}} <: AbstractReductionResult

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

Fields

  • setcovering::SetCovering{Int,T,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.SolutionSizeType
struct SolutionSize{T}

The size of the problem given a configuration.

Fields

  • size: the size of the problem.
  • is_valid: whether the configuration is valid.
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

  • locations::Vector{NTuple{D, T}}: the locations of the vertices
  • radius::Float64: the radius of the unit disk
source
ProblemReductions.cut_sizeMethod
cut_size(g::AbstractGraph, config; weights=UnitWeight(ne(g)))

Return the size of the cut of the graph g with configuration config. The configuration is a vector of boolean numbers as the group indices of vertices. Edges between vertices in different groups are counted as a cut.

source
ProblemReductions.energyMethod
energy(problem::AbstractProblem, config) -> Number

The energy of the problem given the configuration config. Please check the energy_mode for the definition of the energy function.

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.flavor_namesMethod
flavor_names(::Type{<:AbstractProblem}) -> Vector

Returns a vector as the names of the flavors (domain) of a degree of freedom. It falls back to flavors if no method is defined. Use ProblemReductions.name2config and ProblemReductions.config2name to convert between the names and the configuration.

source
ProblemReductions.flavorsMethod
flavors(::Type{<:AbstractProblem}) -> UnitRange
flavors(::GT) where GT<:AbstractProblem -> UnitRange

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

Warning

Flavors is defined a 0:num_flavors-1. To access the previous version of the flavors, use flavor_names.

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_satisfiedMethod
is_satisfied(constraint::LocalConstraint, config) -> Bool
is_satisfied(problem::ConstraintSatisfactionProblem, config) -> Bool

Check if the constraint is satisfied by the configuration config.

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.is_weightedMethod
is_weighted(problem::ConstraintSatisfactionProblem) -> Bool

Check if the problem is weighted. Returns true if the problem has non-unit weights.

source
ProblemReductions.num_flavorsMethod
num_flavors(::Type{<:AbstractProblem}) -> Int
num_flavors(::GT) where GT<:AbstractProblem -> Int

Returns the number of flavors (domain) of a degree of freedom.

source
ProblemReductions.objectivesFunction
objectives(problem::AbstractProblem) -> Vector{<:LocalSolutionSize}

The constraints related to the size of the problem. Each term is associated with weights.

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(problem_type::Type{<:AbstractProblem}, problem::AbstractProblem) -> AbstractReductionResult
reduceto(path::ReductionPath, problem::AbstractProblem) -> ConcatenatedReduction

Reduce the problem problem to a target problem. If the target problem is a single problem type, reduce the problem problem to a target problem of type. Then the result is an instance of AbstractReductionResult. Otherwise, if the target problem is a reduction path, implement a reduction path on a problem. Then the result is a ConcatenatedReduction instance.

Arguments

  • problem_type::Type{<:AbstractProblem} or path::ReductionPath: The target problem type or a reduction path.
  • problem::AbstractProblem: 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_weightsFunction
set_weights(problem::ConstraintSatisfactionProblem, weights) -> ConstraintSatisfactionProblem

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

source
ProblemReductions.solution_sizeMethod
solution_size(problem::AbstractProblem, config) -> SolutionSize

Size of the problem given the configuration config. If you have multiple configurations, use ProblemReductions.solution_size_multiple instead for better performance.

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.variablesMethod
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