Reference

OMEinsumContractionOrders.ExactTreewidthType
const ExactTreewidth{GM} = Treewidth{SafeRules{BT, MMW{3}(), MF}, GM}
ExactTreewidth(; greedy_config = GreedyMethod(nrepeat=1)) = Treewidth(; greedy_config)

ExactTreewidth is a specialization of Treewidth for the SafeRules preprocessing algorithm with the BT elimination algorithm. The BT algorithm is an exact solver for the treewidth problem that implemented in TreeWidthSolver.jl.

source
OMEinsumContractionOrders.GreedyMethodType
GreedyMethod{MT}
GreedyMethod(; α = 0.0, temperature = 0.0, nrepeat=1)

The fast but poor greedy optimizer. Input arguments are

* `α` is the parameter for the loss function, for pairwise interaction, L = size(out) - α * (size(in1) + size(in2))
* `temperature` is the parameter for sampling, if it is zero, the minimum loss is selected; for non-zero, the loss is selected by the Boltzmann distribution, given by p ~ exp(-loss/temperature).
* `nrepeat` is the number of repeatition, returns the best contraction order.
source
OMEinsumContractionOrders.KaHyParBipartiteType
KaHyParBipartite{RT,IT,GM}
KaHyParBipartite(; sc_target, imbalances=collect(0.0:0.005:0.8),
    max_group_size=40, greedy_config=GreedyMethod())

Optimize the einsum code contraction order using the KaHyPar + Greedy approach. This program first recursively cuts the tensors into several groups using KaHyPar, with maximum group size specifed by max_group_size and maximum space complexity specified by sc_target, Then finds the contraction order inside each group with the greedy search algorithm. Other arguments are

  • sc_target is the target space complexity, defined as log2(number of elements in the largest tensor),
  • imbalances is a KaHyPar parameter that controls the group sizes in hierarchical bipartition,
  • max_group_size is the maximum size that allowed to used greedy search,
  • greedy_config is a greedy optimizer.

References

source
OMEinsumContractionOrders.MergeGreedyType
MergeGreedy <: CodeSimplifier
MergeGreedy(; threshhold=-1e-12)

Contraction code simplifier (in order to reduce the time of calling optimizers) that merges tensors greedily if the space complexity of merged tensors is reduced (difference smaller than the threshhold).

source
OMEinsumContractionOrders.SABipartiteType
SABipartite{RT,BT}
SABipartite(; sc_target=25, ntrials=50, βs=0.1:0.2:15.0, niters=1000
    max_group_size=40, greedy_config=GreedyMethod(), initializer=:random)

Optimize the einsum code contraction order using the Simulated Annealing bipartition + Greedy approach. This program first recursively cuts the tensors into several groups using simulated annealing, with maximum group size specifed by max_group_size and maximum space complexity specified by sc_target, Then finds the contraction order inside each group with the greedy search algorithm. Other arguments are

  • size_dict, a dictionary that specifies leg dimensions,
  • sc_target is the target space complexity, defined as log2(number of elements in the largest tensor),
  • max_group_size is the maximum size that allowed to used greedy search,
  • βs is a list of inverse temperature 1/T,
  • niters is the number of iteration in each temperature,
  • ntrials is the number of repetition (with different random seeds),
  • sub_optimizer, the optimizer for the bipartited sub graphs, one can choose GreedyMethod() or TreeSA(),
  • initializer, the partition configuration initializer, one can choose :random or :greedy (slow but better).

References

source
OMEinsumContractionOrders.TreeSAType
TreeSA{RT,IT,GM}
TreeSA(; sc_target=20, βs=collect(0.01:0.05:15), ntrials=10, niters=50,
    sc_weight=1.0, rw_weight=0.2, initializer=:greedy, greedy_config=GreedyMethod(; nrepeat=1))

Optimize the einsum contraction pattern using the simulated annealing on tensor expression tree.

  • sc_target is the target space complexity,
  • ntrials, βs and niters are annealing parameters, doing ntrials indepedent annealings, each has inverse tempteratures specified by βs, in each temperature, do niters updates of the tree.
  • sc_weight is the relative importance factor of space complexity in the loss compared with the time complexity.
  • rw_weight is the relative importance factor of memory read and write in the loss compared with the time complexity.
  • initializer specifies how to determine the initial configuration, it can be :greedy or :random. If it is using :greedy method to generate the initial configuration, it also uses two extra arguments greedy_method and greedy_nrepeat.
  • nslices is the number of sliced legs, default is 0.
  • fixed_slices is a vector of sliced legs, default is [].

References

source
OMEinsumContractionOrders.TreewidthType
struct Treewidth{EL <: EliminationAlgorithm, GM} <: CodeOptimizer
Treewidth(; alg::EL = SafeRules(BT(), MMW{3}(), MF()), greedy_config::GM = GreedyMethod(nrepeat=1))

Tree width based solver. The solvers are implemented in CliqueTrees.jl and TreeWidthSolver.jl. They include:

AlgorithmDescriptionTime ComplexitySpace Complexity
BFSbreadth-first searchO(m + n)O(n)
MCSmaximum cardinality searchO(m + n)O(n)
LexBFSlexicographic breadth-first searchO(m + n)O(m + n)
RCMMDreverse Cuthill-Mckee (minimum degree)O(m + n)O(m + n)
RCMGLreverse Cuthill-Mckee (George-Liu)O(m + n)O(m + n)
MCSMmaximum cardinality search (minimal)O(mn)O(n)
LexMlexicographic breadth-first search (minimal)O(mn)O(n)
AMFapproximate minimum fillO(mn)O(m + n)
MFminimum fillO(mn²)-
MMDmultiple minimum degreeO(mn²)O(m + n)

Detailed descriptions is available in the CliqueTrees.jl.

Fields

  • alg::EL: The algorithm to use for the treewidth calculation. Available elimination algorithms are listed above.
  • greedy_config::GM: The configuration for the greedy method.

Example

julia> optimizer = Treewidth();

julia> eincode = OMEinsumContractionOrders.EinCode([['a', 'b'], ['a', 'c', 'd'], ['b', 'c', 'e', 'f'], ['e'], ['d', 'f']], ['a'])
ab, acd, bcef, e, df -> a

julia> size_dict = Dict([c=>(1<<i) for (i,c) in enumerate(['a', 'b', 'c', 'd', 'e', 'f'])]...)
Dict{Char, Int64} with 6 entries:
  'f' => 64
  'a' => 2
  'c' => 8
  'd' => 16
  'e' => 32
  'b' => 4

julia> optcode = optimize_code(eincode, size_dict, optimizer)
ab, ab -> a
├─ fac, bcf -> ab
│  ├─ df, acd -> fac
│  │  ├─ df
│  │  └─ acd
│  └─ e, bcef -> bcf
│     ├─ e
│     └─ bcef
└─ ab
source
OMEinsumContractionOrders.contraction_complexityMethod
contraction_complexity(eincode, size_dict) -> ContractionComplexity

Returns the time, space and read-write complexity of the einsum contraction. The returned object contains 3 fields:

  • time complexity tc defined as log2(number of element-wise multiplications).
  • space complexity sc defined as log2(size of the maximum intermediate tensor).
  • read-write complexity rwc defined as log2(the number of read-write operations).
source
OMEinsumContractionOrders.flopMethod
flop(eincode, size_dict) -> Int

Returns the number of iterations, which is different with the true floating point operations (FLOP) by a factor of 2.

source
OMEinsumContractionOrders.label_elimination_orderMethod
label_elimination_order(code) -> Vector

Returns a vector of labels sorted by the order they are eliminated in the contraction tree. The contraction tree is specified by code, which e.g. can be a NestedEinsum instance.

source
OMEinsumContractionOrders.optimize_codeFunction
optimize_code(eincode, size_dict, optimizer = GreedyMethod(), simplifier=nothing, permute=true) -> optimized_eincode

Optimize the einsum contraction code and reduce the time/space complexity of tensor network contraction. Returns a NestedEinsum instance. Input arguments are

  • eincode is an einsum contraction code instance, one of DynamicEinCode, StaticEinCode or NestedEinsum.
  • size is a dictionary of "edge label=>edge size" that contains the size information, one can use uniformsize(eincode, 2) to create a uniform size.
  • optimizer is a CodeOptimizer instance, should be one of GreedyMethod, Treewidth, KaHyParBipartite, SABipartite or TreeSA. Check their docstrings for details.
  • simplifier is one of MergeVectors or MergeGreedy.
  • optimize the permutation if permute is true.

Examples

julia> using OMEinsum

julia> code = ein"ij, jk, kl, il->"
ij, jk, kl, il -> 

julia> optimize_code(code, uniformsize(code, 2), TreeSA());
source
OMEinsumContractionOrders.optimize_greedyMethod
optimize_greedy(eincode, size_dict; α = 0.0, temperature = 0.0, nrepeat=1)

Greedy optimizing the contraction order and return a NestedEinsum object. Check the docstring of tree_greedy for detailed explaination of other input arguments.

source
OMEinsumContractionOrders.optimize_kahyparMethod
optimize_kahypar(code, size_dict; sc_target, max_group_size=40, imbalances=0.0:0.01:0.2, greedy_method=MinSpaceOut(), greedy_nrepeat=1)

Optimize the einsum code contraction order using the KaHyPar + Greedy approach. size_dict is a dictionary that specifies leg dimensions. Check the docstring of KaHyParBipartite for detailed explaination of other input arguments.

source
OMEinsumContractionOrders.optimize_kahypar_autoMethod
optimize_kahypar_auto(code, size_dict; max_group_size=40, sub_optimizer = GreedyMethod())

Find the optimal contraction order automatically by determining the sc_target with bisection. It can fail if the tree width of your graph is larger than 100.

source
OMEinsumContractionOrders.optimize_saMethod
optimize_sa(code, size_dict; sc_target, max_group_size=40, βs=0.1:0.2:15.0, niters=1000, ntrials=50,
       sub_optimizer = GreedyMethod(), initializer=:random)

Optimize the einsum code contraction order using the Simulated Annealing bipartition + Greedy approach. size_dict is a dictionary that specifies leg dimensions. Check the docstring of SABipartite for detailed explaination of other input arguments.

References

source
OMEinsumContractionOrders.optimize_treeMethod
optimize_tree(code, size_dict; sc_target=20, βs=0.1:0.1:10, ntrials=2, niters=100, sc_weight=1.0, rw_weight=0.2, initializer=:greedy, greedy_method=MinSpaceOut(), fixed_slices=[])

Optimize the einsum contraction pattern specified by code, and edge sizes specified by size_dict. Check the docstring of TreeSA for detailed explaination of other input arguments.

source
OMEinsumContractionOrders.optimize_treewidthMethod
optimize_treewidth(optimizer, eincode, size_dict)

Optimizing the contraction order via solve the exact tree width of the line graph corresponding to the eincode and return a NestedEinsum object. Check the docstring of treewidth_method for detailed explaination of other input arguments.

source
OMEinsumContractionOrders.tree_greedyMethod
tree_greedy(incidence_list, log2_sizes; α = 0.0, temperature = 0.0, nrepeat=1)

Compute greedy order, and the time and space complexities, the rows of the incidence_list are vertices and columns are edges. log2_sizes are defined on edges. α is the parameter for the loss function, for pairwise interaction, L = size(out) - α * (size(in1) + size(in2)) temperature is the parameter for sampling, if it is zero, the minimum loss is selected; for non-zero, the loss is selected by the Boltzmann distribution, given by p ~ exp(-loss/temperature).

source