Title: | A Toolbox for the Multi-Criteria Minimum Spanning Tree Problem |
---|---|
Description: | Algorithms to approximate the Pareto-front of multi-criteria minimum spanning tree problems. |
Authors: | Jakob Bossek [aut, cre] |
Maintainer: | Jakob Bossek <[email protected]> |
License: | BSD_2_clause + file LICENSE |
Version: | 1.1.1 |
Built: | 2025-02-02 03:11:26 UTC |
Source: | https://github.com/jakobbossek/mcmst |
The mcMST package provides a set of algorithms to approximate the Pareto-optimal set/front of multi-criteria minimum spanning tree (mcMST) problems.
Currently, the following algorithms are included:
A multi-criteria version of Prim's algorithm for the single-objective MST (see [1]).
Evolutionary multi-objective algorithm operating on the Pruefer-encoding as proposed by Zhou and Gen [2].
Evolutionary multi-objective algorithm operating on a direct edge list encoding. This algorithm applies a sub-tree based mutation operator as proposed by Bossek and Grimme [3].
A simple method to enumerate all Pareto-optimal solutions of a given combinatorial problem. This method is not limited to mcMST problems.
[1] Knowles, J. D., and Corne, D. W. 2001. A Comparison of Encodings and Algorithms for Multiobjective Minimum Spanning Tree Problems. In Proceedings of the 2001 Congress on Evolutionary Computation (Ieee Cat. No.01TH8546), 1:544–51 vol. 1. doi:10.1109/CEC.2001.934439.
[2] Zhou, G., and Gen, M. 1999. Genetic Algorithm Approach on Multi-Criteria Minimum Spanning Tree Problem. European Journal of Operational Research 114 (1): 141–52. doi:https://doi.org/10.1016/S0377- 2217(98)00016-2.
[3] Bossek, J., and Grimme, C. 2017. A Pareto-Beneficial Sub-Tree Mutation for the Multi-Criteria Minimum Spanning Tree Problem. In Proceedings of the 2017 IEEE Symposium Series on Computational Intelligence. (accepted)
Convert characteristic vector to edge list.
charVecToEdgelist(charvec)
charVecToEdgelist(charvec)
charvec |
[ |
[matrix
] Edge list.
Other transformation functions:
edgeListToCharVec()
,
nodelistToEdgelist()
,
permutationToCharVec()
,
permutationToEdgelist()
,
prueferToCharVec()
,
prueferToEdgeList()
# here we generate a random Pruefer-code representing # a random spanning tree of a graph with n = 10 nodes pcode = sample(1:10, 8, replace = TRUE)#' edgelist = charVecToEdgelist(prueferToCharVec(pcode))
# here we generate a random Pruefer-code representing # a random spanning tree of a graph with n = 10 nodes pcode = sample(1:10, 8, replace = TRUE)#' edgelist = charVecToEdgelist(prueferToCharVec(pcode))
Given a list of objects and a function which computes a similarity
measure between two objects of the list, computeSimilarity
returns a
similarity matrix.
computeSimilarityMatrix(set, sim.fun, ...)
computeSimilarityMatrix(set, sim.fun, ...)
set |
[ |
sim.fun |
[ |
... |
[any] |
[matrix(n, n)
] matrix with
being the length of
set
.
Convert edge list to characteristic vector.
edgeListToCharVec(edgelist, n = NULL)
edgeListToCharVec(edgelist, n = NULL)
edgelist |
[ |
n |
[ |
[integer
] Characteristic vector cv with cv[i] = 1 if the i-th edge is in the tree.
Other transformation functions:
charVecToEdgelist()
,
nodelistToEdgelist()
,
permutationToCharVec()
,
permutationToEdgelist()
,
prueferToCharVec()
,
prueferToEdgeList()
# first we generate a small edge list by hand # (assume the given graph has n = 4 nodes) edgelist = matrix(c(1, 2, 2, 4, 3, 4), ncol = 3) print(edgelist) # next we transform the edge into # a characteristic vector cvec = edgeListToCharVec(edgelist, n = 4) print(cvec)
# first we generate a small edge list by hand # (assume the given graph has n = 4 nodes) edgelist = matrix(c(1, 2, 2, 4, 3, 4), ncol = 3) print(edgelist) # next we transform the edge into # a characteristic vector cvec = edgeListToCharVec(edgelist, n = 4) print(cvec)
These functions enumerate all candidate solutions for
a certain combinatorial optimization problem, e.g., all permutations
for a TSP or all Pruefer-codes for a MST problem. Note that the output
grows exponentially with the instance size n
.
enumerateTSP(n) enumerateMST(n)
enumerateTSP(n) enumerateMST(n)
n |
[ |
[matrix
] Each row contains a candidate solution.
sols = enumerateTSP(4L) sols = enumerateMST(4L)
sols = enumerateTSP(4L) sols = enumerateMST(4L)
The instance is composed of two
symmetric weight matrices. The first weight is drawn independently at
random from a distribution, the second one
from a
distribution (see references).
genRandomMCGP(n)
genRandomMCGP(n)
n |
[ |
[grapherator
] Graph.
This is a simple wrapper around the much more flexible graph generation system in package grapherator.
Zhou, G. and Gen, M. Genetic Algorithm Approach on Multi-Criteria Minimum Spanning Tree Problem. In: European Journal of Operational Research (1999).
Knowles, JD & Corne, DW 2001, A comparison of encodings and algorithms for multiobjective minimum spanning tree problems. in Proceedings of the IEEE Conference on Evolutionary Computation, ICEC|Proc IEEE Conf Evol Comput Proc ICEC. vol. 1, Institute of Electrical and Electronics Engineers , pp. 544-551, Congress on Evolutionary Computation 2001, Soul, 1 July.
g = genRandomMCGP(10L) ## Not run: pl = grapherator::plot(g) ## End(Not run)
g = genRandomMCGP(10L) ## Not run: pl = grapherator::plot(g) ## End(Not run)
Generate a random spanning tree of a graph given the number of nodes of the problem instance.
genRandomSpanningTree(n, type = "pruefer")
genRandomSpanningTree(n, type = "pruefer")
n |
[ |
type |
[ |
[integer
| matrix(2, n)
] Return type depends on type
.
genRandomSpanningTree(10) genRandomSpanningTree(10, type = "edgelist")
genRandomSpanningTree(10) genRandomSpanningTree(10, type = "edgelist")
Generate a set of random spanning trees of a graph given the number of nodes of the problem instance.
genRandomSpanningTrees(m, n, type = "pruefer", simplify = TRUE)
genRandomSpanningTrees(m, n, type = "pruefer", simplify = TRUE)
m |
[ |
n |
[ |
type |
[ |
simplify |
[ |
[list
| matrix
] Result type depends on simplify
and type
.
genRandomSpanningTrees(3, 10) genRandomSpanningTrees(3, 10, simplify = FALSE) genRandomSpanningTrees(3, 10, type = "edgelist")
genRandomSpanningTrees(3, 10) genRandomSpanningTrees(3, 10, simplify = FALSE) genRandomSpanningTrees(3, 10, type = "edgelist")
Given two spanning trees, the function returns the subtrees of the intersection of these.
getCommonSubtrees(x, y, n = NULL)
getCommonSubtrees(x, y, n = NULL)
x |
[ |
y |
[ |
n |
[ |
[list
] List of matrizes. Each matrix contains the edges of one
connected subtree.
# assume we have a graph with n = 10 nodes n.nodes = 10 # we define two trees (matrices with colwise edges) stree1 = matrix(c(1, 2, 1, 3, 2, 4, 5, 6, 6, 7), byrow = FALSE, nrow = 2) stree2 = matrix(c(1, 3, 1, 2, 2, 4, 5, 8, 6, 7), byrow = FALSE, nrow = 2) # ... and compute all common subtrees subtrees = getCommonSubtrees(stree1, stree2, n = 10)
# assume we have a graph with n = 10 nodes n.nodes = 10 # we define two trees (matrices with colwise edges) stree1 = matrix(c(1, 2, 1, 3, 2, 4, 5, 6, 6, 7), byrow = FALSE, nrow = 2) stree2 = matrix(c(1, 3, 1, 2, 2, 4, 5, 8, 6, 7), byrow = FALSE, nrow = 2) # ... and compute all common subtrees subtrees = getCommonSubtrees(stree1, stree2, n = 10)
Function which expects a problem instance of a combinatorial optimization problem (e.g., MST), a multi-objective function and a solution enumerator, i.e., a function which enumerates all possible solutions (e.g., all Pruefer codes in case of a MST problem) and determines both the Pareto front and Pareto set by exhaustive enumeration.
getExactFront(instance, obj.fun, enumerator.fun, n.objectives, simplify = TRUE)
getExactFront(instance, obj.fun, enumerator.fun, n.objectives, simplify = TRUE)
instance |
[any] |
obj.fun |
[ |
enumerator.fun |
[ |
n.objectives |
[ |
simplify |
[ |
[list
] List with elements pareto.set
(matrix of Pareto-optimal solutions)
and pareto.front
(matrix of corresponding weight vectors).
This method exhaustively enumerates all possible solutions of a given multi-objective combinatorial optimization problem. Thus, it is limited to small input size due to combinatorial explosion.
# here we enumerate all Pareto-optimal solutions of a bi-objective mcMST problem # we use the Pruefer-code enumerator. Thus, we need to define an objective # function, which is able to handle this type of endcoding objfunMCMST = function(pcode, instance) { getWeight(instance, prueferToEdgeList(pcode)) } # next we generate a random bi-objective graph g = genRandomMCGP(5L) # ... and finally compute the exact front of g res = getExactFront(g, obj.fun = objfunMCMST, enumerator.fun = enumerateMST, n.objectives = 2L) ## Not run: plot(res$pareto.front) ## End(Not run)
# here we enumerate all Pareto-optimal solutions of a bi-objective mcMST problem # we use the Pruefer-code enumerator. Thus, we need to define an objective # function, which is able to handle this type of endcoding objfunMCMST = function(pcode, instance) { getWeight(instance, prueferToEdgeList(pcode)) } # next we generate a random bi-objective graph g = genRandomMCGP(5L) # ... and finally compute the exact front of g res = getExactFront(g, obj.fun = objfunMCMST, enumerator.fun = enumerateMST, n.objectives = 2L) ## Not run: plot(res$pareto.front) ## End(Not run)
Internally mcMSTPrim
is called with weights set accordingly.
getExtremeSolutions(graph)
getExtremeSolutions(graph)
graph |
[ |
[matrix(2, 2)
] The i-th column contains the objective vector
of the extreme i-th extreme solution
Makes use of Kirchhoff's matrix tree theorem to compute the number of spanning trees of a given graph in polynomial time.
getNumberOfSpanningTrees(graph)
getNumberOfSpanningTrees(graph)
graph |
[ |
[integer(1)
]
# generate complete graph g = genRandomMCGP(10) # this is equal to 10^8 (Cayley's theorem) getNumberOfSpanningTrees(g)
# generate complete graph g = genRandomMCGP(10) # this is equal to 10^8 (Cayley's theorem) getNumberOfSpanningTrees(g)
Given a grapherator
object this function
returns a random spanning tree. The tree generation process is a simple heuristic:
A random weight from a -distribution is assigned to each edge of the
graph. Next, a spanning tree is computed by
spantree
.
getRandomSpanningTree(graph)
getRandomSpanningTree(graph)
graph |
[ |
[matrix
] Edge list of spanning tree edges.
Most likely this heuristic does not produce each spanning tree with equal probability.
g = genRandomMCGP(10L) stree = getRandomSpanningTree(g)
g = genRandomMCGP(10L) stree = getRandomSpanningTree(g)
Get the overall costs/weight of a subgraph given its edgelist.
getWeight(graph, edgelist, obj.types = NULL)
getWeight(graph, edgelist, obj.types = NULL)
graph |
[ |
edgelist |
[ |
obj.types |
[ |
[numeric(2)
] Weight vector.
# generate a random bi-objective graph g = genRandomMCGP(5) # generate a random Pruefer code, i.e., a random spanning tree of g pcode = sample(1:5, 3, replace = TRUE) getWeight(g, prueferToEdgeList(pcode)) getWeight(g, prueferToEdgeList(pcode), obj.types = "bottleneck")
# generate a random bi-objective graph g = genRandomMCGP(5) # generate a random Pruefer code, i.e., a random spanning tree of g pcode = sample(1:5, 3, replace = TRUE) getWeight(g, prueferToEdgeList(pcode)) getWeight(g, prueferToEdgeList(pcode), obj.types = "bottleneck")
Evolutionary multi-objective algorithm to solve the
multi-objective minimum spanning tree problem. The algorithm relies
to mutation only to generate offspring. The package contains the subgraph mutator
(see mutSubgraphMST
) or a simple one-edge exchange mutator
(see mutEdgeExchange
). Of course, the user may use any
custom mutator which operators on edge lists as well
(see makeMutator
).
mcMSTEmoaBG( instance, mu, lambda = mu, mut = NULL, selMating = NULL, selSurvival = ecr::selNondom, ref.point = NULL, max.iter = 100L, ... )
mcMSTEmoaBG( instance, mu, lambda = mu, mut = NULL, selMating = NULL, selSurvival = ecr::selNondom, ref.point = NULL, max.iter = 100L, ... )
instance |
[ |
mu |
[ |
lambda |
[ |
mut |
[ |
selMating |
[ |
selSurvival |
[ |
ref.point |
[ |
max.iter |
[ |
... |
[ |
[ecr_result
] List of type ecr_result
with the following components:
The ecr_optimization_task
.
Logger object.
Indizes of the non-dominated solutions in the last population.
(n x d) matrix of the approximated non-dominated front where n is the number of non-dominated points and d is the number of objectives.
Matrix of decision space values resulting with objective values given in pareto.front.
Last population.
Character string describing the reason of termination.
Bossek, J., and Grimme, C. A Pareto-Beneficial Sub-Tree Mutation for the Multi-Criteria Minimum Spanning Tree Problem. In Proceedings of the 2017 IEEE Symposium Series on Computational Intelligence (2017). (accepted)
Mutators mutSubgraphMST
and mutEdgeExchange
Other mcMST EMOAs:
mcMSTEmoaZhou()
Other mcMST algorithms:
mcMSTEmoaZhou()
,
mcMSTPrim()
inst = genRandomMCGP(10) res = mcMSTEmoaBG(inst, mu = 20L, max.iter = 100L) print(res$pareto.front) print(tail(getStatistics(res$log)))
inst = genRandomMCGP(10) res = mcMSTEmoaBG(inst, mu = 20L, max.iter = 100L) print(res$pareto.front) print(tail(getStatistics(res$log)))
Evolutionary multi-objective algorithm to solve the
multi-objective minimum spanning tree problem. The algorithm adopts the
so-called Pruefer-number as the encoding for spanning trees. A Pruefer-number
for a graph with nodes is a sequence of
numbers from
. Cayleys theorem states, that a complete graph width n nodes
has exactly
spanning trees.
The algorithm uses mutation only: each component of an individual is replaced
uniformly at random with another node number from the node set.
mcMSTEmoaZhou( instance, mu, lambda = mu, mut = mutUniformPruefer, selMating = ecr::selSimple, selSurvival = ecr::selNondom, ref.point = NULL, max.iter = 100L )
mcMSTEmoaZhou( instance, mu, lambda = mu, mut = mutUniformPruefer, selMating = ecr::selSimple, selSurvival = ecr::selNondom, ref.point = NULL, max.iter = 100L )
instance |
[ |
mu |
[ |
lambda |
[ |
mut |
[ |
selMating |
[ |
selSurvival |
[ |
ref.point |
[ |
max.iter |
[ |
[ecr_result
] List of type ecr_result
with the following components:
The ecr_optimization_task
.
Logger object.
Indizes of the non-dominated solutions in the last population.
(n x d) matrix of the approximated non-dominated front where n is the number of non-dominated points and d is the number of objectives.
Matrix of decision space values resulting with objective values given in pareto.front.
Last population.
Character string describing the reason of termination.
Zhou, G. and Gen, M. Genetic Algorithm Approach on Multi-Criteria Minimum Spanning Tree Problem. In: European Journal of Operational Research (1999).
Mutator mutUniformPruefer
Other mcMST EMOAs:
mcMSTEmoaBG()
Other mcMST algorithms:
mcMSTEmoaBG()
,
mcMSTPrim()
Approximates the Pareto-optimal mcMST front of a multi-objective
graph problem by iteratively applying Prim's algorithm for the single-objective
MST problem to a scalarized version of the problem. I.e., the weight vector
of an edge
is substituted with a weighted
sum
for different weights
.
mcMSTPrim(instance, n.lambdas = NULL, lambdas = NULL)
mcMSTPrim(instance, n.lambdas = NULL, lambdas = NULL)
instance |
[ |
n.lambdas |
[ |
lambdas |
[ |
[list
] List with component pareto.front
.
Note that this procedure can only find socalled supported efficient solutions, i.e., solutions on the convex hull of the Pareto-optimal front.
J. D. Knowles and D. W. Corne, "A comparison of encodings and algorithms for multiobjective minimum spanning tree problems," in Proceedings of the 2001 Congress on Evolutionary Computation (IEEE Cat. No.01TH8546), vol. 1, 2001, pp. 544–551 vol. 1.
Other mcMST algorithms:
mcMSTEmoaBG()
,
mcMSTEmoaZhou()
g = genRandomMCGP(30) res = mcMSTPrim(g, n.lambdas = 50) print(res$pareto.front)
g = genRandomMCGP(30) res = mcMSTPrim(g, n.lambdas = 50) print(res$pareto.front)
Each edge is replaced with another feasible edge with probability p. By default p = 1/m where m is the number of edges, i.e., in expectation one edge is replaced. The operators maintains the spanning tree property, i.e., the resulting edge list is indeed the edge list of a spanning tree.
mutEdgeExchange(ind, p = 1/ncol(ind), instance = NULL)
mutEdgeExchange(ind, p = 1/ncol(ind), instance = NULL)
ind |
[ |
p |
[ |
instance |
[ |
[matrix(2, m)
] Mutated edge list.
Evolutionary multi-objective algorithm mcMSTEmoaBG
Other mcMST EMOA mutators:
mutKEdgeExchange()
,
mutSubforestMST()
,
mutSubgraphMST()
,
mutUniformPruefer()
Let be the number of spanning tree edges. Then, the operator
selects
edges randomly and replaces each of the
edges with another feasible edge.
mutKEdgeExchange(ind, k = 1L, instance = NULL)
mutKEdgeExchange(ind, k = 1L, instance = NULL)
ind |
[ |
k |
[ |
instance |
[ |
Evolutionary multi-objective algorithm mcMSTEmoaBG
Other mcMST EMOA mutators:
mutEdgeExchange()
,
mutSubforestMST()
,
mutSubgraphMST()
,
mutUniformPruefer()
mutForestMST
drops k edges randomly. In consequence the
tree is decomposed into k+1 subtrees (forest). Now the operator reconnects the
subtrees by constructing a minimum spanning tree between the components.
mutSubforestMST(ind, sigma = ncol(ind), scalarize = FALSE, instance = NULL)
mutSubforestMST(ind, sigma = ncol(ind), scalarize = FALSE, instance = NULL)
ind |
[ |
sigma |
[ |
scalarize |
[ |
instance |
[ |
[matrix(2, m)
] Mutated edge list.
Evolutionary multi-objective algorithm mcMSTEmoaBG
Other mcMST EMOA mutators:
mutEdgeExchange()
,
mutKEdgeExchange()
,
mutSubgraphMST()
,
mutUniformPruefer()
mutSubgraphMST
selects a random edge e = (u, v) and traverses
the tree starting form u and v respectively until a connected subtree of at most
sigma
edges is selected. Then the subtree is replaced with the optimal spanning subtree
regarding one of the objectives with equal probability.
mutSubgraphMST( ind, sigma = floor(ncol(ind)/2), scalarize = FALSE, instance = NULL )
mutSubgraphMST( ind, sigma = floor(ncol(ind)/2), scalarize = FALSE, instance = NULL )
ind |
[ |
sigma |
[ |
scalarize |
[ |
instance |
[ |
[matrix(2, m)
] Mutated edge list.
Evolutionary multi-objective algorithm mcMSTEmoaBG
Other mcMST EMOA mutators:
mutEdgeExchange()
,
mutKEdgeExchange()
,
mutSubforestMST()
,
mutUniformPruefer()
mutUniformPruefer
replaces each component of a Pruefer code of length n - 2
with probability p
with a random node number between 1 and n.
mutUniformPruefer(ind, p = 1/length(ind))
mutUniformPruefer(ind, p = 1/length(ind))
ind |
[ |
p |
[ |
[integer
] Mutated Pruefer code.
Evolutionary multi-objective algorithm mcMSTEmoaZhou
Other mcMST EMOA mutators:
mutEdgeExchange()
,
mutKEdgeExchange()
,
mutSubforestMST()
,
mutSubgraphMST()
Convert sequence of nodes to edge list.
nodelistToEdgelist(nodelist)
nodelistToEdgelist(nodelist)
nodelist |
[ |
[matrix
] Edge list.
Other transformation functions:
charVecToEdgelist()
,
edgeListToCharVec()
,
permutationToCharVec()
,
permutationToEdgelist()
,
prueferToCharVec()
,
prueferToEdgeList()
# first generate a random permutation, e.g., representing # a roundtrip tour in a graph nodelist = sample(1:8) # now convert into an edge list nodelistToEdgelist(nodelist)
# first generate a random permutation, e.g., representing # a roundtrip tour in a graph nodelist = sample(1:8) # now convert into an edge list nodelistToEdgelist(nodelist)
Convert permutation to characteristic vector.
permutationToCharVec(perm, n)
permutationToCharVec(perm, n)
perm |
[ |
n |
[ |
[integer
] Characteristic vector cv with cv[i] = 1 if the i-th edge is in the tree.
Other transformation functions:
charVecToEdgelist()
,
edgeListToCharVec()
,
nodelistToEdgelist()
,
permutationToEdgelist()
,
prueferToCharVec()
,
prueferToEdgeList()
# first generate a random permutation, e.g., representing # a roundtrip tour in a graph perm = sample(1:10) print(perm) # now convert into an edge list permutationToCharVec(perm, n = 10)
# first generate a random permutation, e.g., representing # a roundtrip tour in a graph perm = sample(1:10) print(perm) # now convert into an edge list permutationToCharVec(perm, n = 10)
Convert permutation to edge list.
permutationToEdgelist(perm)
permutationToEdgelist(perm)
perm |
[ |
[matrix(2, length(perm))
] Edge list.
Other transformation functions:
charVecToEdgelist()
,
edgeListToCharVec()
,
nodelistToEdgelist()
,
permutationToCharVec()
,
prueferToCharVec()
,
prueferToEdgeList()
# first generate a random permutation, e.g., representing # a roundtrip tour in a graph perm = sample(1:10) print(perm) # now convert into an edge list permutationToEdgelist(perm)
# first generate a random permutation, e.g., representing # a roundtrip tour in a graph perm = sample(1:10) print(perm) # now convert into an edge list permutationToEdgelist(perm)
Given a list of graphs and a list of solutions (encoded as edge lists) for each graph the function generates each one plot. This is a 2D-scatterplot of edge weights of the graph. Size and colour of each point indicate the number of solutions the edge is part of.
plotEdgeFrequency(graphs, approx.sets, facet.args = list(), names = NULL)
plotEdgeFrequency(graphs, approx.sets, facet.args = list(), names = NULL)
graphs |
[ |
approx.sets |
[ |
facet.args |
[ |
names |
[ |
[ggplot
]
Other result visualization:
plotEdges()
g = genRandomMCGP(50L) res = mcMSTEmoaBG(mu = 10L, max.iter = 50, instance = g, scalarize = TRUE) ## Not run: plotEdgeFrequency(list(g), list(res$pareto.set)) ## End(Not run)
g = genRandomMCGP(50L) res = mcMSTEmoaBG(mu = 10L, max.iter = 50, instance = g, scalarize = TRUE) ## Not run: plotEdgeFrequency(list(g), list(res$pareto.set)) ## End(Not run)
Given a list of characteristic vectors (graphs) the function plots an embedding of the nodes in the Euclidean plane and depicts an edge if and only if it is contained in at least one of the graphs. The edge thickness indicates the number of graphs the edge is part of.
plotEdges(x, n = NULL, normalize = TRUE, ...)
plotEdges(x, n = NULL, normalize = TRUE, ...)
x |
[ |
n |
[ |
normalize |
[ |
... |
[any] |
Nothing
Other result visualization:
plotEdgeFrequency()
Convert Pruefer code to characteristic vector.
prueferToCharVec(pcode)
prueferToCharVec(pcode)
pcode |
[ |
[integer
] Characteristic vector cv with cv[i] = 1 if the i-th edge is in the tree.
Other transformation functions:
charVecToEdgelist()
,
edgeListToCharVec()
,
nodelistToEdgelist()
,
permutationToCharVec()
,
permutationToEdgelist()
,
prueferToEdgeList()
# here we generate a random Pruefer-code representing # a random spanning tree of a graph with n = 10 nodes pcode = sample(1:10, 8, replace = TRUE) print(pcode) print(prueferToCharVec(pcode))
# here we generate a random Pruefer-code representing # a random spanning tree of a graph with n = 10 nodes pcode = sample(1:10, 8, replace = TRUE) print(pcode) print(prueferToCharVec(pcode))
Convert Pruefer code to edge list.
prueferToEdgeList(pcode)
prueferToEdgeList(pcode)
pcode |
[ |
[matrix(2, length(pcode) + 1)
] Edge list.
Other transformation functions:
charVecToEdgelist()
,
edgeListToCharVec()
,
nodelistToEdgelist()
,
permutationToCharVec()
,
permutationToEdgelist()
,
prueferToCharVec()
# here we generate a random Pruefer-code representing # a random spanning tree of a graph with n = 10 nodes pcode = sample(1:10, 8, replace = TRUE) print(pcode) edgelist = prueferToEdgeList(pcode) print(edgelist)
# here we generate a random Pruefer-code representing # a random spanning tree of a graph with n = 10 nodes pcode = sample(1:10, 8, replace = TRUE) print(pcode) edgelist = prueferToEdgeList(pcode) print(edgelist)
Sample random weights
for weighted-sum scalarization.
sampleWeights(n)
sampleWeights(n)
n |
[ |
[numeric
] Weight vector.
sampleWeights(2) weights = replicate(10, sampleWeights(3L)) colSums(weights)
sampleWeights(2) weights = replicate(10, sampleWeights(3L)) colSums(weights)
Given a list of weight matrizes weight.mats
and a vector
of numeric weights, the function returns a single weight matrix. Each component
of the resulting matrix is the weighted sum of the corresponding components of
the weight matrizes passed.
scalarizeWeights(weight.mats, lambdas)
scalarizeWeights(weight.mats, lambdas)
weight.mats |
[ |
lambdas |
[ |
[matrix
]
Functions which expect two (spanning) trees and return a measure
of similiarity between those. Function getNumberOfCommonEdges
returns
the (normalized) number of shared edges and function getSizeOfLargestCommonSubtree
returns the (normalized) size of the largest connected subtree which is located in
both trees.
getNumberOfCommonEdges(x, y, n = NULL, normalize = TRUE) getSizeOfLargestCommonSubtree(x, y, n = NULL, normalize = TRUE)
getNumberOfCommonEdges(x, y, n = NULL, normalize = TRUE) getSizeOfLargestCommonSubtree(x, y, n = NULL, normalize = TRUE)
x |
[ |
y |
[ |
n |
[ |
normalize |
[ |
[numeric(1)
] Measure
# Here we generate two random spanning trees of a complete # graph with 10 nodes set.seed(1) st1 = prueferToEdgeList(sample(1:10, size = 8, replace = TRUE)) st2 = prueferToEdgeList(sample(1:10, size = 8, replace = TRUE)) # Now check the number of common edges NCE = getNumberOfCommonEdges(st1, st2) # And the size of the largest common subtree SLS = getSizeOfLargestCommonSubtree(st1, st2)
# Here we generate two random spanning trees of a complete # graph with 10 nodes set.seed(1) st1 = prueferToEdgeList(sample(1:10, size = 8, replace = TRUE)) st2 = prueferToEdgeList(sample(1:10, size = 8, replace = TRUE)) # Now check the number of common edges NCE = getNumberOfCommonEdges(st1, st2) # And the size of the largest common subtree SLS = getSizeOfLargestCommonSubtree(st1, st2)