Random graphs¶
The methods defined here appear in sage.graphs.graph_generators
.
- sage.graphs.generators.random.RandomBarabasiAlbert(n, m, seed=None)¶
Return a random graph created using the Barabasi-Albert preferential attachment model.
A graph with \(m\) vertices and no edges is initialized, and a graph of \(n\) vertices is grown by attaching new vertices each with \(m\) edges that are attached to existing vertices, preferentially with high degree.
INPUT:
n
– number of vertices in the graphm
– number of edges to attach from each new nodeseed
– arandom.Random
seed or a Pythonint
for the random number generator (default:None
)
EXAMPLES:
We show the edge list of a random graph on 6 nodes with \(m = 2\):
sage: G = graphs.RandomBarabasiAlbert(6,2) sage: G.order(), G.size() (6, 8) sage: G.degree_sequence() # random [4, 3, 3, 2, 2, 2]
We plot a random graph on 12 nodes with \(m = 3\):
sage: ba = graphs.RandomBarabasiAlbert(12,3) sage: ba.show() # long time
We view many random graphs using a graphics array:
sage: g = [] sage: j = [] sage: for i in range(1,10): ....: k = graphs.RandomBarabasiAlbert(i+3, 3) ....: g.append(k) sage: for i in range(3): ....: n = [] ....: for m in range(3): ....: n.append(g[3*i + m].plot(vertex_size=50, vertex_labels=False)) ....: j.append(n) sage: G = graphics_array(j) sage: G.show() # long time
When \(m = 1\), the generated graph is a tree:
sage: graphs.RandomBarabasiAlbert(6, 1).is_tree() True
- sage.graphs.generators.random.RandomBicubicPlanar(n)¶
Return the graph of a random bipartite cubic map with \(3 n\) edges.
INPUT:
\(n\) – an integer (at least \(1\))
OUTPUT:
a graph with multiple edges (no embedding is provided)
The algorithm used is described in [Sch1999]. This samples a random rooted bipartite cubic map, chosen uniformly at random.
First one creates a random binary tree with \(n\) vertices. Next one turns this into a blossoming tree (at random) and reads the contour word of this blossoming tree.
Then one performs a rotation on this word so that this becomes a balanced word. There are three ways to do that, one is picked at random. Then a graph is build from the balanced word by iterated closure (adding edges).
In the returned graph, the three edges incident to any given vertex are colored by the integers 0, 1 and 2.
See also
the auxiliary method
blossoming_contour()
EXAMPLES:
sage: n = randint(200, 300) sage: G = graphs.RandomBicubicPlanar(n) sage: G.order() == 2*n True sage: G.size() == 3*n True sage: G.is_bipartite() and G.is_planar() and G.is_regular(3) True sage: dic = {'red':[v for v in G.vertices() if v[0] == 'n'], ....: 'blue': [v for v in G.vertices() if v[0] != 'n']} sage: G.plot(vertex_labels=False,vertex_size=20,vertex_colors=dic) Graphics object consisting of ... graphics primitives
- sage.graphs.generators.random.RandomBipartite(n1, n2, p, set_position=False)¶
Returns a bipartite graph with \(n1+n2\) vertices such that any edge from \([n1]\) to \([n2]\) exists with probability \(p\).
INPUT:
n1, n2
– Cardinalities of the two setsp
– Probability for an edge to existset_position
– boolean (defaultFalse
); if set toTrue
, we assign positions to the vertices so that the set of cardinality \(n1\) is on the line \(y=1\) and the set of cardinality \(n2\) is on the line \(y=0\).
EXAMPLES:
sage: g = graphs.RandomBipartite(5, 2, 0.5) sage: g.vertices() [(0, 0), (0, 1), (0, 2), (0, 3), (0, 4), (1, 0), (1, 1)]
- sage.graphs.generators.random.RandomBlockGraph(m, k, kmax=None, incidence_structure=False)¶
Return a Random Block Graph.
A block graph is a connected graph in which every biconnected component (block) is a clique.
See also
Wikipedia article Block_graph for more details on these graphs
is_block_graph()
– test if a graph is a block graph
INPUT:
m
– integer; number of blocks (at least one).k
– integer; minimum number of vertices of a block (at least two).kmax
– integer (default:None
) By default, each block has \(k\) vertices. When the parameter \(kmax\) is specified (with \(kmax \geq k\)), the number of vertices of each block is randomly chosen between \(k\) and \(kmax\).incidence_structure
– boolean (default:False
) when set toTrue
, the incidence structure of the graphs is returned instead of the graph itself, that is the list of the lists of vertices in each block. This is useful for the creation of some hypergraphs.
OUTPUT:
A Graph when
incidence_structure==False
(default), and otherwise an incidence structure.EXAMPLES:
A block graph with a single block is a clique:
sage: B = graphs.RandomBlockGraph(1, 4) sage: B.is_clique() True
A block graph with blocks of order 2 is a tree:
sage: B = graphs.RandomBlockGraph(10, 2) sage: B.is_tree() True
Every biconnected component of a block graph is a clique:
sage: B = graphs.RandomBlockGraph(5, 3, kmax=6) sage: blocks,cuts = B.blocks_and_cut_vertices() sage: all(B.is_clique(block) for block in blocks) True
A block graph with blocks of order \(k\) has \(m*(k-1)+1\) vertices:
sage: m, k = 6, 4 sage: B = graphs.RandomBlockGraph(m, k) sage: B.order() == m*(k-1)+1 True
Test recognition methods:
sage: B = graphs.RandomBlockGraph(6, 2, kmax=6) sage: B.is_block_graph() True sage: B in graph_classes.Block True
Asking for the incidence structure:
sage: m, k = 6, 4 sage: IS = graphs.RandomBlockGraph(m, k, incidence_structure=True) sage: from sage.combinat.designs.incidence_structures import IncidenceStructure sage: IncidenceStructure(IS) Incidence structure with 19 points and 6 blocks sage: m*(k-1)+1 19
- sage.graphs.generators.random.RandomBoundedToleranceGraph(n)¶
Returns a random bounded tolerance graph.
The random tolerance graph is built from a random bounded tolerance representation by using the function \(ToleranceGraph\). This representation is a list \(((l_0,r_0,t_0), (l_1,r_1,t_1), ..., (l_k,r_k,t_k))\) where \(k = n-1\) and \(I_i = (l_i,r_i)\) denotes a random interval and \(t_i\) a random positive value less then or equal to the length of the interval \(I_i\). The width of the representation is limited to n**2 * 2**n.
Note
The tolerance representation used to create the graph can be recovered using
get_vertex()
orget_vertices()
.INPUT:
n
– number of vertices of the random graph.
EXAMPLES:
Every (bounded) tolerance graph is perfect. Hence, the chromatic number is equal to the clique number
sage: g = graphs.RandomBoundedToleranceGraph(8) sage: g.clique_number() == g.chromatic_number() True
- sage.graphs.generators.random.RandomChordalGraph(n, algorithm='growing', k=None, l=None, f=None, s=None)¶
Return a random chordal graph of order
n
.A Graph \(G\) is said to be chordal if it contains no induced hole (a cycle of length at least 4). Equivalently, \(G\) is chordal if it has a perfect elimination orderings, if each minimal separator is a clique, or if it is the intersection graphs of subtrees of a tree. See the Wikipedia article Chordal_graph.
This generator implements the algorithms proposed in [SHET2018] for generating random chordal graphs as the intersection graph of \(n\) subtrees of a tree of order \(n\).
The returned graph is not necessarily connected.
INPUT:
n
– integer; the number of nodes of the graphalgorithm
– string (default:"growing"
); the choice of the algorithm for randomly selecting \(n\) subtrees of a random tree of order \(n\). Possible choices are:"growing"
– for each subtree \(T_i\), the algorithm picks a size \(k_i\) randomly from \([1,k]\). Then a random node of \(T\) is chosen as the first node of \(T_i\). In each of the subsequent \(k_i - 1\) iterations, it picks a random node in the neighborhood of \(T_i\) and adds it to \(T_i\)."connecting"
– for each subtree \(T_i\), it first selects \(k_i\) nodes of \(T\), where \(k_i\) is a random integer from a Poisson distribution with mean \(l\). \(T_i\) is then generated to be the minimal subtree containing the selected \(k_i\) nodes. This implies that a subtree will most likely have many more nodes than those selected initially, and this must be taken into consideration when choosing \(l\)."pruned"
– for each subtree \(T_i\), it randomly selects a fraction \(f\) of the edges on the tree and removes them. The number of edges to delete, say \(l\), is calculated as \(\lfloor (n - 1) f \rfloor\), which will leave \(l + 1\) subtrees in total. Then, it determines the sizes of the \(l + 1\) subtrees and stores the distinct values. Finally, it picks a random size \(k_i\) from the set of largest \(100(1-s)\%\) of distinct values, and randomly chooses a subtree with size \(k_i\).
k
– integer (default:None
); maximum size of a subtree. If not specified (None
), the maximum size is set to \(\sqrt{n}\). This parameter is used only whenalgorithm="growing"
. Seegrowing_subtrees()
for more details.l
– a strictly positive real number (default:None
); mean of a Poisson distribution. If not specified, the mean in set to \(\log_2{n}\). This parameter is used only whenalgorithm="connecting"
. Seeconnecting_nodes()
for more details.f
– a rational number (default:None
); the edge deletion fraction. This value must be chosen in \([0..1]\). If not specified, this parameter is set to \(\frac{1}{n-1}\). This parameter is used only whenalgorithm="pruned"
. Seepruned_tree()
for more details.s
– a real number between 0 and 1 (default:None
); selection barrier for the size of trees. If not specified, this parameter is set to \(0.5\). This parameter is used only whenalgorithm="pruned"
. Seepruned_tree()
for more details.
EXAMPLES:
sage: from sage.graphs.generators.random import RandomChordalGraph sage: T = RandomChordalGraph(20, algorithm="growing", k=5) sage: T.is_chordal() True sage: T = RandomChordalGraph(20, algorithm="connecting", l=3) sage: T.is_chordal() True sage: T = RandomChordalGraph(20, algorithm="pruned", f=1/3, s=.5) sage: T.is_chordal() True
See also
growing_subtrees()
connecting_nodes()
pruned_tree()
- sage.graphs.generators.random.RandomGNM(n, m, dense=False, seed=None)¶
Returns a graph randomly picked out of all graphs on n vertices with m edges.
INPUT:
n
- number of vertices.m
- number of edges.dense
- whether to use NetworkX’s dense_gnm_random_graph or gnm_random_graphseed
- arandom.Random
seed or a Pythonint
for the random number generator (default:None
).
EXAMPLES: We show the edge list of a random graph on 5 nodes with 10 edges.
sage: graphs.RandomGNM(5, 10).edges(labels=False) [(0, 1), (0, 2), (0, 3), (0, 4), (1, 2), (1, 3), (1, 4), (2, 3), (2, 4), (3, 4)]
We plot a random graph on 12 nodes with m = 12.
sage: gnm = graphs.RandomGNM(12, 12) sage: gnm.show() # long time
We view many random graphs using a graphics array:
sage: g = [] sage: j = [] sage: for i in range(9): ....: k = graphs.RandomGNM(i+3, i^2-i) ....: g.append(k) sage: for i in range(3): ....: n = [] ....: for m in range(3): ....: n.append(g[3*i + m].plot(vertex_size=50, vertex_labels=False)) ....: j.append(n) sage: G = graphics_array(j) sage: G.show() # long time
- sage.graphs.generators.random.RandomGNP(n, p, seed=None, fast=True, algorithm='Sage')¶
Returns a random graph on \(n\) nodes. Each edge is inserted independently with probability \(p\).
INPUT:
n
– number of nodes of the graphp
– probability of an edgeseed
- arandom.Random
seed or a Pythonint
for the random number generator (default:None
).fast
– boolean set to True (default) to use the algorithm with time complexity in \(O(n+m)\) proposed in [BB2005a]. It is designed for generating large sparse graphs. It is faster than other algorithms for LARGE instances (try it to know whether it is useful for you).algorithm
– By default (`algorithm='Sage'
), this function uses the algorithm implemented in`sage.graphs.graph_generators_pyx.pyx
. Whenalgorithm='networkx'
, this function calls the NetworkX functionfast_gnp_random_graph
, unlessfast=False
, thengnp_random_graph
. Try them to know which algorithm is the best for you. Thefast
parameter is not taken into account by the ‘Sage’ algorithm so far.
REFERENCES:
PLOTTING: When plotting, this graph will use the default spring-layout algorithm, unless a position dictionary is specified.
EXAMPLES: We show the edge list of a random graph on 6 nodes with probability \(p = .4\):
sage: set_random_seed(0) sage: graphs.RandomGNP(6, .4).edges(labels=False) [(0, 1), (0, 5), (1, 2), (2, 4), (3, 4), (3, 5), (4, 5)]
We plot a random graph on 12 nodes with probability \(p = .71\):
sage: gnp = graphs.RandomGNP(12,.71) sage: gnp.show() # long time
We view many random graphs using a graphics array:
sage: g = [] sage: j = [] sage: for i in range(9): ....: k = graphs.RandomGNP(i+3,.43) ....: g.append(k) sage: for i in range(3): ....: n = [] ....: for m in range(3): ....: n.append(g[3*i + m].plot(vertex_size=50, vertex_labels=False)) ....: j.append(n) sage: G = graphics_array(j) sage: G.show() # long time sage: graphs.RandomGNP(4,1) Complete graph: Graph on 4 vertices
- sage.graphs.generators.random.RandomHolmeKim(n, m, p, seed=None)¶
Return a random graph generated by the Holme and Kim algorithm for graphs with power law degree distribution and approximate average clustering.
INPUT:
n
– number of verticesm
– number of random edges to add for each new nodep
– probability of adding a triangle after adding a random edgeseed
– arandom.Random
seed or a Pythonint
for the random number generator (default:None
)
From the NetworkX documentation: the average clustering has a hard time getting above a certain cutoff that depends on \(m\). This cutoff is often quite low. Note that the transitivity (fraction of triangles to possible triangles) seems to go down with network size. It is essentially the Barabasi-Albert growth model with an extra step that each random edge is followed by a chance of making an edge to one of its neighbors too (and thus a triangle). This algorithm improves on B-A in the sense that it enables a higher average clustering to be attained if desired. It seems possible to have a disconnected graph with this algorithm since the initial \(m\) nodes may not be all linked to a new node on the first iteration like the BA model.
EXAMPLES:
We check that a random graph on 8 nodes with 2 random edges per node and a probability \(p = 0.5\) of forming triangles contains a triangle:
sage: G = graphs.RandomHolmeKim(8, 2, 0.5) sage: G.order(), G.size() (8, 12) sage: C3 = graphs.CycleGraph(3) sage: G.subgraph_search(C3) Subgraph of (): Graph on 3 vertices
sage: G = graphs.RandomHolmeKim(12, 3, .3) sage: G.show() # long time
REFERENCE:
- sage.graphs.generators.random.RandomIntervalGraph(n)¶
Returns a random interval graph.
An interval graph is built from a list \((a_i,b_i)_{1\leq i \leq n}\) of intervals : to each interval of the list is associated one vertex, two vertices being adjacent if the two corresponding intervals intersect.
A random interval graph of order \(n\) is generated by picking random values for the \((a_i,b_j)\), each of the two coordinates being generated from the uniform distribution on the interval \([0,1]\).
This definitions follows [BF2001].
Note
The vertices are named 0, 1, 2, and so on. The intervals used to create the graph are saved with the graph and can be recovered using
get_vertex()
orget_vertices()
.INPUT:
n
(integer) – the number of vertices in the random graph.
EXAMPLES:
As for any interval graph, the chromatic number is equal to the clique number
sage: g = graphs.RandomIntervalGraph(8) sage: g.clique_number() == g.chromatic_number() True
- sage.graphs.generators.random.RandomLobster(n, p, q, seed=None)¶
Returns a random lobster.
A lobster is a tree that reduces to a caterpillar when pruning all leaf vertices. A caterpillar is a tree that reduces to a path when pruning all leaf vertices (q=0).
INPUT:
n
- expected number of vertices in the backbonep
- probability of adding an edge to the backboneq
- probability of adding an edge (claw) to the armsseed
- arandom.Random
seed or a Pythonint
for the random number generator (default:None
).
EXAMPLES:
We check a random graph with 12 backbone nodes and probabilities \(p = 0.7\) and \(q = 0.3\):
sage: G = graphs.RandomLobster(12, 0.7, 0.3) sage: leaves = [v for v in G.vertices() if G.degree(v) == 1] sage: G.delete_vertices(leaves) # caterpillar sage: leaves = [v for v in G.vertices() if G.degree(v) == 1] sage: G.delete_vertices(leaves) # path sage: s = G.degree_sequence() sage: if G: ....: assert s[-2:] == [1, 1] ....: assert all(d == 2 for d in s[:-2])
sage: G = graphs.RandomLobster(9, .6, .3) sage: G.show() # long time
- sage.graphs.generators.random.RandomNewmanWattsStrogatz(n, k, p, seed=None)¶
Return a Newman-Watts-Strogatz small world random graph on \(n\) vertices.
From the NetworkX documentation: first create a ring over \(n\) nodes. Then each node in the ring is connected with its \(k\) nearest neighbors. Then shortcuts are created by adding new edges as follows: for each edge \(u-v\) in the underlying “\(n\)-ring with \(k\) nearest neighbors”; with probability \(p\) add a new edge \(u-w\) with randomly-chosen existing node \(w\). In contrast with
networkx.watts_strogatz_graph()
, no edges are removed.INPUT:
n
– number of verticesk
– each vertex is connected to its \(k\) nearest neighborsp
– the probability of adding a new edge for each edgeseed
– arandom.Random
seed or a Pythonint
for the random number generator (default:None
)
EXAMPLES:
We check that the generated graph contains a cycle of order \(n\):
sage: G = graphs.RandomNewmanWattsStrogatz(7, 2, 0.2) sage: G.order() 7 sage: C7 = graphs.CycleGraph(7) sage: G.subgraph_search(C7) Subgraph of (): Graph on 7 vertices sage: G.diameter() <= C7.diameter() True
sage: G = graphs.RandomNewmanWattsStrogatz(12, 2, .3) sage: G.show() # long time
REFERENCE:
- sage.graphs.generators.random.RandomRegular(d, n, seed=None)¶
Return a random \(d\)-regular graph on \(n\) vertices, or
False
on failure.Since every edge is incident to two vertices, \(n\times d\) must be even.
INPUT:
d
– degreen
– number of verticesseed
– arandom.Random
seed or a Pythonint
for the random number generator (default:None
)
EXAMPLES:
We check that a random graph with 8 nodes each of degree 3 is 3-regular:
sage: G = graphs.RandomRegular(3, 8) sage: G.is_regular(k=3) True sage: G.degree_histogram() [0, 0, 0, 8]
sage: G = graphs.RandomRegular(3, 20) sage: if G: ....: G.show() # random output, long time
REFERENCES:
- sage.graphs.generators.random.RandomRegularBipartite(n1, n2, d1, set_position=False)¶
Return a random regular bipartite graph on \(n1 + n2\) vertices.
The bipartite graph has \(n1 * d1\) edges. Hence, \(n2\) must divide \(n1 * d1\). Each vertex of the set of cardinality \(n1\) has degree \(d1\) (which can be at most \(n2\)) and each vertex in the set of cardinality \(n2\) has degree \((n1 * d1) / n2\). The bipartite graph has no multiple edges.
This generator implements an algorithm inspired by that of [MW1990] for the uniform generation of random regular bipartite graphs. It performs well when \(d1 = o(n2^{1/3})\) or (\(n2 - d1 = o(n2^{1/3})\)). In other cases, the running time can be huge. Note that the currently implemented algorithm does not generate uniformly random graphs.
INPUT:
n1, n2
– number of vertices in each sided1
– degree of the vertices in the set of cardinality \(n1\).set_position
– boolean (defaultFalse
); if set toTrue
, we assign positions to the vertices so that the set of cardinality \(n1\) is on the line \(y=1\) and the set of cardinality \(n2\) is on the line \(y=0\).
EXAMPLES:
sage: g = graphs.RandomRegularBipartite(4, 6, 3) sage: g.order(), g.size() (10, 12) sage: set(g.degree()) {2, 3} sage: graphs.RandomRegularBipartite(1, 2, 2, set_position=True).get_pos() {0: (1, 1.0), 1: (0, 0), 2: (2.0, 0.0)} sage: graphs.RandomRegularBipartite(2, 1, 1, set_position=True).get_pos() {0: (0, 1), 1: (2.0, 1.0), 2: (1, 0.0)} sage: graphs.RandomRegularBipartite(2, 3, 3, set_position=True).get_pos() {0: (0, 1), 1: (3.0, 1.0), 2: (0, 0), 3: (1.5, 0.0), 4: (3.0, 0.0)} sage: graphs.RandomRegularBipartite(2, 3, 3, set_position=False).get_pos()
- sage.graphs.generators.random.RandomShell(constructor, seed=None)¶
Return a random shell graph for the constructor given.
INPUT:
constructor
– a list of 3-tuples \((n, m, d)\), each representing a shell, where:n
– the number of vertices in the shellm
– the number of edges in the shelld
– the ratio of inter (next) shell edges to intra shell edges
seed
– arandom.Random
seed or a Pythonint
for the random number generator (default:None
)
EXAMPLES:
sage: G = graphs.RandomShell([(10,20,0.8),(20,40,0.8)]) sage: G.order(), G.size() (30, 52) sage: G.show() # long time
- sage.graphs.generators.random.RandomToleranceGraph(n)¶
Returns a random tolerance graph.
The random tolerance graph is built from a random tolerance representation by using the function \(ToleranceGraph\). This representation is a list \(((l_0,r_0,t_0), (l_1,r_1,t_1), ..., (l_k,r_k,t_k))\) where \(k = n-1\) and \(I_i = (l_i,r_i)\) denotes a random interval and \(t_i\) a random positive value. The width of the representation is limited to n**2 * 2**n.
Note
The vertices are named 0, 1, …, n-1. The tolerance representation used to create the graph is saved with the graph and can be recovered using
get_vertex()
orget_vertices()
.INPUT:
n
– number of vertices of the random graph.
EXAMPLES:
Every tolerance graph is perfect. Hence, the chromatic number is equal to the clique number
sage: g = graphs.RandomToleranceGraph(8) sage: g.clique_number() == g.chromatic_number() True
- sage.graphs.generators.random.RandomTree(n)¶
Returns a random tree on \(n\) nodes numbered \(0\) through \(n-1\).
By Cayley’s theorem, there are \(n^{n-2}\) trees with vertex set \(\{0,1,...,n-1\}\). This constructor chooses one of these uniformly at random.
ALGORITHM:
The algorithm works by generating an \((n-2)\)-long random sequence of numbers chosen independently and uniformly from \(\{0,1,\ldots,n-1\}\) and then applies an inverse Prufer transformation.
INPUT:
n
- number of vertices in the tree
EXAMPLES:
sage: G = graphs.RandomTree(10) sage: G.is_tree() True sage: G.show() # long time
- sage.graphs.generators.random.RandomTreePowerlaw(n, gamma=3, tries=1000, seed=None)¶
Return a tree with a power law degree distribution, or
False
on failure.From the NetworkX documentation: a trial power law degree sequence is chosen and then elements are swapped with new elements from a power law distribution until the sequence makes a tree (size = order - 1).
INPUT:
n
– number of verticesgamma
– exponent of power law distributiontries
– number of attempts to adjust sequence to make a treeseed
– arandom.Random
seed or a Pythonint
for the random number generator (default:None
)
EXAMPLES:
We check that the generated graph is a tree:
sage: G = graphs.RandomTreePowerlaw(10, 3) sage: G.is_tree() True sage: G.order(), G.size() (10, 9)
sage: G = graphs.RandomTreePowerlaw(15, 2) sage: if G: ....: G.show() # random output, long time
- sage.graphs.generators.random.RandomTriangulation(n, set_position=False, k=3)¶
Return a random inner triangulation of an outer face of degree
k
withn
vertices in total.An inner triangulation is a plane graph all of whose faces (except the outer/unbounded face) are triangles (3-cycles).
INPUT:
n
– the number of vertices of the graphk
– the size of the outer faceset_position
– boolean (defaultFalse
); if set toTrue
, this will compute coordinates for a planar drawing of the graph.
OUTPUT:
A random graph chosen uniformly among the inner triangulations of a rooted \(k\)-gon with \(n\) vertices (including the \(k\) vertices from the outer face). This is a planar graph and comes with a combinatorial embedding. The vertices of the root edge are labelled
-1
and-2
and the outer face is the face returned byGraph.faces()
in which-1
and-2
are consecutive vertices in this order.Because some triangulations have nontrivial automorphism groups, this may not be equal to the uniform distribution among inner triangulations of unrooted \(k\)-gons.
ALGORITHM:
The algorithm is taken from [PS2006], Section 5.
Starting from a planar \(k\)-gonal forest (represented by its contour as a sequence of vertices), one performs local closures, until no one is possible. A local closure amounts to replace in the cyclic contour word a sequence
in1, in2, in3, lf, in3
byin1, in3
.At every step of the algorithm, newly created edges are recorded in a graph, which will be returned at the end. The combinatorial embedding is also computed and recorded in the output graph.
See also
EXAMPLES:
sage: G = graphs.RandomTriangulation(6, True); G Graph on 6 vertices sage: G.is_planar() True sage: G.girth() 3 sage: G.plot(vertex_size=0, vertex_labels=False) Graphics object consisting of 13 graphics primitives sage: H = graphs.RandomTriangulation(7, k=5) sage: sorted(len(f) for f in H.faces()) [3, 3, 3, 3, 3, 3, 3, 5]
- sage.graphs.generators.random.blossoming_contour(t, shift=0)¶
Return a random blossoming of a binary tree \(t\), as a contour word.
This is doing several things simultaneously:
complete the binary tree, by adding leaves labelled
xb
,add a vertex labelled
n
at the middle of every inner edge, with a leaf labelledx
either on the left or on the right (at random),number all vertices (but not leaves) by integers starting from \(shift\),
compute the counter-clockwise contour word of the result.
Initial vertices receive the label
i
.This is an auxiliary function, used for the generation of random planar bicubic maps.
INPUT:
\(t\) – a binary tree (non-empty)
shift
– an integer (default \(0\)), used as a starting index
OUTPUT:
contour word of a random blossoming of \(t\)
EXAMPLES:
sage: from sage.graphs.generators.random import blossoming_contour sage: print(blossoming_contour(BinaryTrees(1).an_element())) [('i', 0), ('xb',), ('i', 0), ('xb',), ('i', 0)] sage: t = BinaryTrees(2).random_element() sage: print(blossoming_contour(t)) # random [('i', 0), ('xb',), ('i', 0), ('n', 2), ('i', 1), ('xb',), ('i', 1), ('xb',), ('i', 1), ('n', 2), ('x',), ('n', 2), ('i', 0)] sage: w = blossoming_contour(BinaryTrees(3).random_element()); len(w) 21 sage: w.count(('xb',)) 4 sage: w.count(('x',)) 2
- sage.graphs.generators.random.connecting_nodes(T, l)¶
Return a list of the vertex sets of
n
randomly chosen subtrees ofT
.This method is part of
RandomChordalGraph()
.ALGORITHM:
For each subtree \(T_i\), we first select \(k_i\) nodes of \(T\), where \(k_i\) is a random integer from a Poisson distribution with mean \(l\). \(T_i\) is then generated to be the minimal subtree that contains the selected \(k_i\) nodes. This implies that a subtree will most likely have many more nodes than those selected initially, and this must be taken into consideration when choosing \(l\).
See [SHET2018] for more details.
INPUT:
T
– a treel
– a strictly positive real number; mean of a Poisson distribution
EXAMPLES:
sage: from sage.graphs.generators.random import connecting_nodes sage: T = graphs.RandomTree(10) sage: S = connecting_nodes(T, 5) sage: len(S) 10
- sage.graphs.generators.random.growing_subtrees(T, k)¶
Return a list of the vertex sets of
n
randomly chosen subtrees ofT
.For a tree of order \(n\), the collection contains \(n\) subtrees with maximum order \(k\) and average order \(\frac{k + 1}{2}\).
This method is part of
RandomChordalGraph()
.ALGORITHM:
For each subtree \(T_i\), the algorithm picks a size \(k_i\) randomly from \([1,k]\). Then a random node of \(T\) is chosen as the first node of \(T_i\). In each of the subsequent \(k_i - 1\) iterations, it picks a random node in the neighborhood of \(T_i\) and adds it to \(T_i\).
See [SHET2018] for more details.
INPUT:
T
– a treek
– a strictly positive integer; maximum size of a subtree
EXAMPLES:
sage: from sage.graphs.generators.random import growing_subtrees sage: T = graphs.RandomTree(10) sage: S = growing_subtrees(T, 5) sage: len(S) 10
- sage.graphs.generators.random.pruned_tree(T, f, s)¶
Return a list of the vertex sets of
n
randomly chosen subtrees ofT
.This method is part of
RandomChordalGraph()
.ALGORITHM:
For each subtree \(T_i\), it randomly selects a fraction \(f\) of the edges on the tree and removes them. The number of edges to delete, say \(l\), is calculated as \(\lfloor((n - 1)f\rfloor\), which will leave \(l + 1\) subtrees in total. Then, it determines the sizes of the \(l + 1\) subtrees and stores the distinct values. Finally, it picks a random size \(k_i\) from the set of largest \(100(1-s)\%\) of distinct values, and randomly chooses a subtree with size \(k_i\).
See [SHET2018] for more details.
INPUT:
T
– a treef
– a rational number; the edge deletion fraction. This value must be chosen in \([0..1]\).s
– a real number between 0 and 1; selection barrier for the size of trees
EXAMPLES:
sage: from sage.graphs.generators.random import pruned_tree sage: T = graphs.RandomTree(11) sage: S = pruned_tree(T, 1/10, 0.5) sage: len(S) 11