Library Reference
Lower Bound Algorithms
CliqueTrees.LowerBoundAlgorithm
— TypeLowerBoundAlgorithm
An algorithm for computing a lower bound to the treewidth of a graph. The options are
type | name | time | space |
---|---|---|---|
MMW | minor-min-width |
CliqueTrees.MMW
— TypeMMW <: LowerBoundAlgorithm
MMW()
The minor-min-width heuristic.
References
- Gogate, Vibhav, and Rina Dechter. "A complete anytime algorithm for treewidth." Proceedings of the 20th conference on Uncertainty in artificial intelligence. 2004.
- Bodlaender, Hans, Thomas Wolle, and Arie Koster. "Contraction and treewidth lower bounds." Journal of Graph Algorithms and Applications 10.1 (2006): 5-49.
CliqueTrees.lowerbound
— Functionlowerbound([weights, ]graph;
alg::WidthOrAlgorithm=DEFAULT_LOWER_BOUND_ALGORITHM)
Compute a lower bound to the treewidth of a graph.
CliqueTrees.DEFAULT_LOWER_BOUND_ALGORITHM
— ConstantDEFAULT_LOWER_BOUND_ALGORITHM = MMW()
Elimination Algorithms
CliqueTrees.EliminationAlgorithm
— TypeEliminationAlgorithm
A graph elimination algorithm computes a permutation of the vertices of a graph, which induces a chordal completion of the graph. The algorithms below generally seek to minimize the fill (number of edges) or width (largest clique) of the completed graph.
type | name | time | space |
---|---|---|---|
BFS | breadth-first search | O(m + n) | O(n) |
MCS | maximum cardinality search | O(m + n) | O(n) |
LexBFS | lexicographic breadth-first search | O(m + n) | O(m + n) |
RCMMD | reverse Cuthill-Mckee (minimum degree) | O(m + n) | O(m + n) |
RCMGL | reverse Cuthill-Mckee (George-Liu) | O(m + n) | O(m + n) |
MCSM | maximum cardinality search (minimal) | O(mn) | O(n) |
LexM | lexicographic breadth-first search (minimal) | O(mn) | O(n) |
AMF | approximate minimum fill | O(mn) | O(m + n) |
MF | minimum fill | O(mn²) | |
MMD | multiple minimum degree | O(mn²) | O(m + n) |
MinimalChordal | MinimalChordal | ||
CompositeRotations | elimination tree rotation | O(m + n) | O(m + n) |
SafeRules | treewith-safe rule-based reduction | ||
SafeSeparators | treewith-safe separator decomposition | ||
ConnectedComponents | connected component decomposition |
The following additional algorithms are implemented as package extensions and require loading an additional package.
type | name | time | space | package |
---|---|---|---|---|
AMD | approximate minimum degree | O(mn) | O(m + n) | AMD.jl |
SymAMD | column approximate minimum degree | O(mn) | O(m + n) | AMD.jl |
METIS | multilevel nested dissection | Metis.jl | ||
Spectral | spectral ordering | Laplacians.jl | ||
FlowCutter | FlowCutter | FlowCutterPACE17_jll.jl | ||
BT | Bouchitte-Todinca | TreeWidthSolver.jl | ||
SAT{libpicosat_jll} | SAT encoding (picosat) | libpicosat_jll.jl | ||
SAT{PicoSAT_jll} | SAT encoding (picosat) | PicoSAT_jll.jl | ||
SAT{Lingeling_jll} | SAT encoding (lingeling) | Lingeling_jll.jl | ||
SAT{CryptoMiniSat_jll} | SAT encoding (cryptominisat) | CryptoMiniSat_jll.jl |
Triangulation Recognition Heuristics
These algorithms are guaranteed to compute perfect elimination orderings for chordal graphs. MCSM
and LexM
are variants of MCS
and LexBFS
that compute minimal orderings. The Lex algorithms were pubished first, and the MCS algorithms were introducd later as simplications. In practice, these algorithms work poorly on non-chordal graphs.
Bandwidth and Envelope Reduction Heuristics
These algorithms seek to minimize the bandwidth and profile of the permuted graph, quantities that upper bound the width and fill of theinduced chordal completion. RCMMD
and RCMGL
are two variants of the reverse Cuthill-McKee algorithm, a type of breadth-first search. They differ in in their choice of starting vertex. In practice, these algorithms work better than the triangulation recognition heuristics and worse than the greedy heuristics.
Greedy Heuristics
These algorithms simulate the elimination process, greedity selecting vertices to eliminate. MMD
selects a vertex of minimum degree, and MF
selects a vertex that induces the least fill. Updating the degree or fill of every vertex after elimination is costly; the algorithms AMD
, SymAMD
, and AMF
are relaxations that work by approximating these values. The AMD
algorithm is the state-of-the-practice for sparse matrix ordering.
Exact Treewidth Algorithms
These algorithm minimizes the treewidth of the completed graph.
This is an NP-hard problem. I recommend wrapping exact treewidth algorithms with preprocessors like SafeRules
or ConnectedComponents
.
Meta Algorithms
These algorithms are parametrized by another algorithm and work by transforming its input or output.
CliqueTrees.PermutationOrAlgorithm
— TypePermutationOrAlgorithm = Union{AbstractVector, EliminationAlgorithm}
Either a permutation or an algorithm.
CliqueTrees.DEFAULT_ELIMINATION_ALGORITHM
— ConstantDEFAULT_ELIMINATION_ALGORITHM = MMD()
The default algorithm.
CliqueTrees.BFS
— TypeBFS <: EliminationAlgorithm
BFS()
The breadth-first search algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = BFS()
BFS
julia> treewidth(graph; alg)
2
CliqueTrees.MCS
— TypeMCS <: EliminationAlgorithm
MCS()
The maximum cardinality search algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = MCS()
MCS
julia> treewidth(graph; alg)
3
References
- Tarjan, Robert E., and Mihalis Yannakakis. "Simple linear-time algorithms to test chordality of graphs, test acyclicity of hypergraphs, and selectively reduce acyclic hypergraphs." SIAM Journal on Computing 13.3 (1984): 566-579.
CliqueTrees.LexBFS
— TypeLexBFS <: EliminationAlgorithm
LexBFS()
The lexicographic breadth-first-search algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = LexBFS()
LexBFS
julia> treewidth(graph; alg)
2
References
- Rose, Donald J., R. Endre Tarjan, and George S. Lueker. "Algorithmic aspects of vertex elimination on graphs." SIAM Journal on Computing 5.2 (1976): 266-283.
CliqueTrees.RCMMD
— TypeRCMMD{A} <: EliminationAlgorithm
RCMMD(alg::Algorithm)
RCMMD()
The reverse Cuthill-McKee algorithm. An initial vertex is selected using the minimum degree heuristic.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = RCMMD(QuickSort)
RCMMD:
Base.Sort.QuickSortAlg()
julia> treewidth(graph; alg)
3
Parameters
alg
: sorting algorithm
References
- Cuthill, Elizabeth, and James McKee. "Reducing the bandwidth of sparse symmetric matrices." Proceedings of the 1969 24th National Conference. 1969.
CliqueTrees.RCMGL
— TypeRCMGL{A} <: EliminationAlgorithm
RCMGL(alg::Algorithm)
RCMGL()
The reverse Cuthill-McKee algorithm. An initial vertex is selected using George and Liu's variant of the GPS algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = RCMGL(QuickSort)
RCMGL:
Base.Sort.QuickSortAlg()
julia> treewidth(graph; alg)
3
Parameters
alg
: sorting algorithm
References
- Cuthill, Elizabeth, and James McKee. "Reducing the bandwidth of sparse symmetric matrices." Proceedings of the 1969 24th National Conference. 1969.
- George, Alan, and Joseph WH Liu. "An implementation of a pseudoperipheral node finder." ACM Transactions on Mathematical Software (TOMS) 5.3 (1979): 284-295.
CliqueTrees.RCM
— TypeRCM = RCMGL
The default variant of the reverse Cuthill-Mckee algorithm.
CliqueTrees.LexM
— TypeLexM <: EliminationAlgorithm
LexM()
A minimal variant of the lexicographic breadth-first-search algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = LexM()
LexM
julia> treewidth(graph; alg)
2
References
- Rose, Donald J., R. Endre Tarjan, and George S. Lueker. "Algorithmic aspects of vertex elimination on graphs." SIAM Journal on Computing 5.2 (1976): 266-283.
CliqueTrees.MCSM
— TypeMCSM <: EliminationAlgorithm
MCSM()
A minimal variant of the maximal cardinality search algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = MCSM()
MCSM
julia> treewidth(graph; alg)
2
References
- Berry, Anne, et al. "Maximum cardinality search for computing minimal triangulations of graphs." Algorithmica 39 (2004): 287-298.
CliqueTrees.AMD
— TypeAMD <: EliminationAlgorithm
AMD(; dense=10.0, aggressive=1.0)
The approximate minimum degree algorithm.
julia> using CliqueTrees
julia> import AMD as AMDLib
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = AMD(; dense=5.0, aggressive=2.0)
AMD:
dense: 5.0
aggressive: 2.0
julia> treewidth(graph; alg)
2
Parameters
dense
: dense row parameteraggressive
: aggressive absorption
References
- Amestoy, Patrick R., Timothy A. Davis, and Iain S. Duff. "An approximate minimum degree ordering algorithm." SIAM Journal on Matrix Analysis and Applications 17.4 (1996): 886-905.
CliqueTrees.SymAMD
— TypeSymAMD <: EliminationAlgorithm
SymAMD(; dense_row=10.0, dense_col=10.0, aggressive=1.0)
The column approximate minimum degree algorithm.
julia> using CliqueTrees
julia> import AMD as AMDLib
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = SymAMD(; dense_row=5.0, dense_col=5.0, aggressive=2.0)
SymAMD:
dense_row: 5.0
dense_col: 5.0
aggressive: 2.0
julia> treewidth(graph, alg)
2
Parameters
dense_row
: dense row parameterdense_column
: dense column parameteraggressive
: aggressive absorption
References
- Davis, Timothy A., et al. "A column approximate minimum degree ordering algorithm." ACM Transactions on Mathematical Software (TOMS) 30.3 (2004): 353-376.
CliqueTrees.AMF
— TypeAMF <: EliminationAlgorithm
AMF(; speed=1)
The approximate minimum fill algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = AMF(; speed=2)
AMF:
speed: 2
julia> treewidth(graph; alg)
2
Parameters
speed
: fill approximation strategy (1
,2
, or,3
)
References
- Rothberg, Edward, and Stanley C. Eisenstat. "Node selection strategies for bottom-up sparse matrix ordering." SIAM Journal on Matrix Analysis and Applications 19.3 (1998): 682-695.
CliqueTrees.MF
— TypeMF <: EliminationAlgorithm
MF()
The greedy minimum fill algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = MF
MF
julia> treewidth(graph; alg)
2
References
- Tinney, William F., and John W. Walker. "Direct solutions of sparse network equations by optimally ordered triangular factorization." Proceedings of the IEEE 55.11 (1967): 1801-1809.
CliqueTrees.MMD
— TypeMMD <: EliminationAlgorithm
MMD(; delta=0)
The multiple minimum degree algorithm.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = MMD(; delta=1)
MMD:
delta: 1
julia> treewidth(graph; alg)
2
Parameters
delta
: tolerance for multiple elimination
References
- Liu, Joseph WH. "Modification of the minimum-degree algorithm by multiple elimination." ACM Transactions on Mathematical Software (TOMS) 11.2 (1985): 141-153.
CliqueTrees.METIS
— TypeMETIS <: EliminationAlgorithm
METIS(; ctype=-1, rtype=-1, nseps=-1, niter=-1, seed=-1,
compress=-1, ccorder=-1, pfactor=-1, ufactor=-1)
The multilevel nested dissection algorithm implemented in METIS.
julia> using CliqueTrees, Metis
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = METIS(ctype=Metis.METIS_CTYPE_RM)
METIS:
ctype: 0
rtype: -1
nseps: -1
niter: -1
seed: -1
compress: -1
ccorder: -1
pfactor: -1
ufactor: -1
julia> treewidth(graph; alg)
3
Parameters
ctype
: matching scheme to be used during coarseningrtype
: algorithm used for refinementnseps
: number of different separators computed at each level of nested dissectionniter
: number of iterations for refinement algorithm at each stage of the uncoarsening processseed
: random seedcompress
: whether to combine vertices with identical adjacency listsccorder
: whether to order connected components separatelypfactor
: minimum degree of vertices that will be ordered lastufactor
: maximum allowed load imbalance partitions
References
- Karypis, George, and Vipin Kumar. "A fast and high quality multilevel scheme for partitioning irregular graphs." SIAM Journal on Scientific Computing 20.1 (1998): 359-392.
CliqueTrees.Spectral
— TypeSpectral <: EliminationAlgorithm
Spectral(; tol=0.0)
The spectral ordering algorithm only works on connected graphs. In order to use it, import the package Laplacians.
julia> using CliqueTrees, Laplacians
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = Spectral(; tol=0.001)
Spectral:
tol: 0.001
julia> treewidth(graph; alg)
4
Parameters
tol
: tolerance for convergence
References
- Barnard, Stephen T., Alex Pothen, and Horst D. Simon. "A spectral algorithm for envelope reduction of sparse matrices." Proceedings of the 1993 ACM/IEEE Conference on Supercomputing. 1993.
CliqueTrees.FlowCutter
— TypeFlowCutter <: EliminationAlgorithm
FlowCutter(; time=5, seed=0)
The FlowCutter algorithm.
julia> using CliqueTrees, FlowCutterPACE17_jll
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = FlowCutter(; time=2, seed=1)
FlowCutter:
time: 2
seed: 1
julia> treewidth(graph; alg)
2
Parameters
time
: run timeseed
: random seed
References
- Strasser, Ben. "Computing tree decompositions with flowcutter: PACE 2017 submission." arXiv preprint arXiv:1709.08949 (2017).
CliqueTrees.BT
— TypeBT <: EliminationAlgorithm
BT()
The Bouchitte-Todinca algorithm.
julia> using CliqueTrees, TreeWidthSolver
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = BT()
BT
julia> treewidth(graph; alg)
2
References
- Korhonen, Tuukka, Jeremias Berg, and Matti Järvisalo. "Solving Graph Problems via Potential Maximal Cliques: An Experimental Evaluation of the Bouchitté-Todinca Algorithm." Journal of Experimental Algorithmics (JEA) 24 (2019): 1-19.
CliqueTrees.SAT
— TypeSAT{H, L, U} <: EliminationAlgorithm
SAT{H}(lb::WidthOrAlgorithn, ub::PermutationOrAlgorithm)
SAT{H}()
Compute a minimum-treewidth permutation using a SAT solver.
julia> using CliqueTrees, libpicosat_jll, PicoSAT_jll, CryptoMiniSat_jll, Lingeling_jll
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = SAT{libpicosat_jll}(MMW(), MF()) # picosat
SAT{libpicosat_jll, MMW, MF}:
MMW
MF
julia> alg = SAT{PicoSAT_jll}(MMW(), MF()) # picosat
SAT{PicoSAT_jll, MMW, MF}:
MMW
MF
julia> alg = SAT{CryptoMiniSat_jll}(MMW(), MF()) # cryptominisat
SAT{CryptoMiniSat_jll, MMW, MF}:
MMW
MF
julia> alg = SAT{Lingeling_jll}(MMW(), MF()) # lingeling
SAT{Lingeling_jll, MMW, MF}:
MMW
MF
julia> treewidth(graph; alg)
2
Parameters
lb
: lower bound algorithmub
: upper bound algorithm
References
- Samer, Marko, and Helmut Veith. "Encoding treewidth into SAT." Theory and Applications of Satisfiability Testing-SAT 2009: 12th International Conference, SAT 2009, Swansea, UK, June 30-July 3, 2009. Proceedings 12. Springer Berlin Heidelberg, 2009.
- Berg, Jeremias, and Matti Järvisalo. "SAT-based approaches to treewidth computation: An evaluation." 2014 IEEE 26th international conference on tools with artificial intelligence. IEEE, 2014.
- Bannach, Max, Sebastian Berndt, and Thorsten Ehlers. "Jdrasil: A modular library for computing tree decompositions." 16th International Symposium on Experimental Algorithms (SEA 2017). Schloss Dagstuhl–Leibniz-Zentrum fuer Informatik, 2017.
CliqueTrees.MinimalChordal
— TypeMinimalChordal{A} <: EliminationAlgorithm
MinimalChordal(alg::PermutationOrAlgorithm)
MinimalChordal()
Evaluate an elimination algorithm, and them improve its output using the MinimalChordal algorithm. The result is guaranteed to be minimal.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg1 = MCS()
MCS
julia> alg2 = MinimalChordal(MCS())
MinimalChordal{MCS}:
MCS
julia> label1, tree1 = cliquetree(graph; alg=alg1);
julia> label2, tree2 = cliquetree(graph; alg=alg2);
julia> FilledGraph(tree1) # more edges
{8, 12} FilledGraph{Int64, Int64}
julia> FilledGraph(tree2) # fewer edges
{8, 11} FilledGraph{Int64, Int64}
Parameters
alg
: elimination algorithm
References
- Blair, Jean RS, Pinar Heggernes, and Jan Arne Telle. "A practical algorithm for making filled graphs minimal." Theoretical Computer Science 250.1-2 (2001): 125-141.
CliqueTrees.CompositeRotations
— TypeCompositeRotations{C, A} <: EliminationAlgorithm
CompositeRotations(clique::AbstractVector, alg::EliminationAlgorithm)
CompositeRotations(clique::AbstractVector)
Evaluate an eliminaton algorithm, ensuring that the given clique is at the end of the ordering.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg = CompositeRotations([2], MCS())
CompositeRotations{Vector{Int64}, MCS}:
clique: [2]
MCS
julia> order, index = permutation(graph; alg);
julia> order # 2 is the last vertex in the ordering
8-element Vector{Int64}:
4
5
7
8
3
6
1
2
Parameters
clique
: cliquealg
: elimination algorithm
References
- Liu, Joseph WH. "Equivalent sparse matrix reordering by elimination tree rotations." Siam Journal on Scientific and Statistical Computing 9.3 (1988): 424-444.
CliqueTrees.SafeRules
— TypeSafeRules{L, U} <: EliminationAlgorithm
SafeRules(lb::WidthOrAlgorithm, ub::EliminationAlgororithm)
SafeRules()
Preprocess a graph using safe reduction rules. The algorithm lb
is used to compute a lower bound to the treewidth; better lower bounds allow the algorithm to perform more reductions.
julia> using CliqueTrees, TreeWidthSolver
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> alg1 = BT()
BT
julia> alg2 = SafeRules(MMW(), BT())
SafeRules{MMW, BT}:
MMW
BT
julia> @time treewidth(graph; alg=alg1) # slower
0.000177 seconds (1.41 k allocations: 90.031 KiB)
2
julia> @time treewidth(graph; alg=alg2) # faster
0.000066 seconds (237 allocations: 16.875 KiB)
2
Parameters
lb
: lower bound algorithmub
: elimination algorithm
References
- Bodlaender, Hans L., et al. "Pre-processing for triangulation of probabilistic networks." (2001).
- Bodlaender, Hans L., Arie M.C.A. Koster, and Frank van den Eijkhof. "Preprocessing rules for triangulation of probabilistic networks." Computational Intelligence 21.3 (2005): 286-305.
- van den Eijkhof, Frank, Hans L. Bodlaender, and Arie M.C.A. Koster. "Safe reduction rules for weighted treewidth." Algorithmica 47 (2007): 139-158.
CliqueTrees.SafeSeparators
— TypeSafeSeparators{M, A} <: EliminationAlgorithm
SafeSeparators(min::PermutationOrAlgorithm, alg::EliminationAlgorithm)
SeparatorReducton()
Apple an elimination algorithm to the atoms of an almost-clique separator decomposition. The algorithm min
is used to compute the decomposition.
The algorithm min
must compute a minimimal ordering. I recommend wrapping a heuristic algorithm with MinimalChordal
.
Parameters
min
: minimal elimination algorithmalg
: elimination algorithm
References
- Bodlaender, Hans L., and Arie MCA Koster. "Safe separators for treewidth." Discrete Mathematics 306.3 (2006): 337-350.
- Tamaki, Hisao. "A heuristic for listing almost-clique minimal separators of a graph." arXiv preprint arXiv:2108.07551 (2021).
CliqueTrees.ConnectedComponents
— TypeConnectedComponents{A} <: EliminationAlgorithm
ConnectedComponents(alg::PermutationOrAlgorithm)
ConnectedComponents()
Apply an elimination algorithm to each connected component of a graph.
Parameters
alg
: elimination algorithm
CliqueTrees.permutation
— Functionpermutation([weights, ]graph;
alg::PermutationOrAlgorithm=DEFAULT_ELIMINATION_ALGORITHM)
Construct a fill-reducing permutation of the vertices of a simple graph.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> order, index = permutation(graph);
julia> order
8-element Vector{Int64}:
4
1
2
8
5
3
6
7
julia> index == invperm(order)
true
CliqueTrees.mcs
— Functionmcs(graph[, clique])
Perform a maximum cardinality search, optionally specifying a clique to be ordered last. Returns the inverse permutation.
Supernodes
CliqueTrees.SupernodeType
— TypeSupernodeType
A type of supernode partition. The options are
type | name |
---|---|
Nodal | nodal supernode partition |
Maximal | maximal supernode partition |
Fundamental | fundamental supernode partition |
CliqueTrees.DEFAULT_SUPERNODE_TYPE
— ConstantDEFAULT_SUPERNODE_TYPE = Maximal()
The default supernode partition.
CliqueTrees.Nodal
— TypeNodal <: SupernodeType
A nodal supernode partition.
CliqueTrees.Maximal
— TypeMaximal <: SupernodeType
A maximal supernode partition.
CliqueTrees.Fundamental
— TypeFundamental <: SupernodeType
A fundamental supernode partition.
Linked Lists
CliqueTrees.SinglyLinkedList
— TypeSinglyLinkedList{I, Head, Next} <: AbstractLinkedList{I}
SinglyLinkedList{I}(n::Integer)
A singly linked list of distinct natural numbers. This type supports the iteration interface.
julia> using CliqueTrees
julia> list = SinglyLinkedList{Int}(10)
SinglyLinkedList{Int64, Array{Int64, 0}, Vector{Int64}}:
julia> pushfirst!(list, 4, 5, 6, 7, 8, 9)
SinglyLinkedList{Int64, Array{Int64, 0}, Vector{Int64}}:
4
5
6
7
8
⋮
julia> collect(list)
6-element Vector{Int64}:
4
5
6
7
8
9
Trees
CliqueTrees.AbstractTree
— TypeAbstractTree{V} = Union{Tree{V}, SupernodeTree{V}, CliqueTree{V}}
A rooted forest. This type implements the indexed tree interface.
CliqueTrees.rootindices
— Functionrootindices(tree::AbstractTree)
Get the roots of a rooted forest.
CliqueTrees.firstchildindex
— Functionfirstchildindex(tree::AbstractTree, i::Integer)
Get the first child of node i
. Returns nothing
if i
is a leaf.
CliqueTrees.ancestorindices
— Functionancestorindices(tree::AbstractTree, i::Integer)
Get the proper ancestors of node i
.
Trees
CliqueTrees.Tree
— TypeTree{V} <: AbstractUnitRange{V}
Tree(tree::AbstractTree)
Tree{V}(tree::AbstractTree) where V
A rooted forest with vertices of type V
. This type implements the indexed tree interface.
CliqueTrees.eliminationtree
— Functioneliminationtree([weights, ]graph;
alg::PermutationOrAlgorithm=DEFAULT_ELIMINATION_ALGORITHM)
Construct a tree-depth decomposition of a simple graph.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> label, tree = eliminationtree(graph);
julia> tree
8-element Tree{Int64}:
8
└─ 7
├─ 5
└─ 6
├─ 1
├─ 3
│ └─ 2
└─ 4
Supernode Trees
CliqueTrees.SupernodeTree
— TypeSupernodeTree{V} <: AbstractVector{UnitRange{V}}
A supernodal elimination tree with vertices of type V
. This type implements the indexed tree interface.
CliqueTrees.supernodetree
— Functionsupernodetree(graph;
alg::PermutationOrAlgorithm=DEFAULT_ELIMINATION_ALGORITHM,
snd::SupernodeType=DEFAULT_SUPERNODE_TYPE)
Construct a supernodal elimination tree.
Clique Trees
CliqueTrees.Clique
— TypeClique{V, E} <: AbstractVector{V}
A clique of a clique tree.
CliqueTrees.CliqueTree
— TypeCliqueTree{V, E} <: AbstractVector{Clique{V, E}}
A clique tree with vertices of type V
and edges of type E
. This type implements the indexed tree interface.
CliqueTrees.cliquetree
— Functioncliquetree([weights, ]graph;
alg::PermutationOrAlgorithm=DEFAULT_ELIMINATION_ALGORITHM,
snd::SupernodeType=DEFAULT_SUPERNODE_TYPE)
Construct a tree decomposition of a simple graph. The vertices of the graph are first ordered by a fill-reducing permutation computed by the algorithm alg
. The size of the resulting decomposition is determined by the supernode partition snd
.
julia> using CliqueTrees
julia> graph = [
0 1 0 0 0 0 0 0
1 0 1 0 0 1 0 0
0 1 0 1 0 1 1 1
0 0 1 0 0 0 0 0
0 0 0 0 0 1 1 0
0 1 1 0 1 0 0 0
0 0 1 0 1 0 0 1
0 0 1 0 0 0 1 0
];
julia> label, tree = cliquetree(graph);
julia> tree
6-element CliqueTree{Int64, Int64}:
[6, 7, 8]
└─ [5, 7, 8]
├─ [1, 5]
├─ [3, 5, 7]
│ └─ [2, 3]
└─ [4, 5, 8]
CliqueTrees.treewidth
— Functiontreewidth([weights, ]tree::CliqueTree)
Compute the width of a clique tree.
treewidth([weights, ]graph;
alg::PermutationOrAlgorithm=DEFAULT_ELIMINATION_ALGORITHM)
Compute an upper bound to the tree width of a simple graph.
CliqueTrees.separator
— Functionseparator(clique::Clique)
Get the separator of a clique.
separator(tree::CliqueTree, i::Integer)
Get the separator at node i
.
CliqueTrees.residual
— Functionresidual(clique::Clique)
Get the residual of a clique.
residual(tree::CliqueTree, i::Integer)
Get the residual at node i
.
Filled Graphs
CliqueTrees.FilledGraph
— TypeFilledGraph{V, E} <: AbstractGraph{V}
A filled graph.
CliqueTrees.ischordal
— Functionischordal(graph)
Determine whether a simple graph is chordal.
CliqueTrees.isperfect
— Functionisperfect(graph, order::AbstractVector[, index::AbstractVector])
Determine whether an fill-reducing permutation is perfect.