1% This file is part of the Stanford GraphBase (c) Stanford University 1993
2@i boilerplate.w %<< legal stuff: PLEASE READ IT BEFORE MAKING ANY CHANGES!
3@i gb_types.w
4
5\def\title{GB\_\,BASIC}
6
7\prerequisite{GB\_\,GRAPH}
8@* Introduction. This GraphBase module contains six subroutines that generate
9standard graphs of various types, together with six routines that combine or
10transform existing graphs.
11
12Simple examples of the use of these routines can be found in the
13demonstration programs {\sc QUEEN} and {\sc QUEEN\_WRAP}.
14
15@<gb_basic.h@>=
16extern Graph *board(); /* moves on generalized chessboards */
17extern Graph *simplex(); /* generalized triangular configurations */
18extern Graph *subsets(); /* patterns of subset intersection */
19extern Graph *perms(); /* permutations of a multiset */
20extern Graph *parts(); /* partitions of an integer */
21extern Graph *binary(); /* binary trees */
22@#
23extern Graph *complement(); /* the complement of a graph */
24extern Graph *gunion(); /* the union of two graphs */
25extern Graph *intersection(); /* the intersection of two graphs */
26extern Graph *lines(); /* the line graph of a graph */
27extern Graph *product(); /* the product of two graphs */
28extern Graph *induced(); /* a graph induced from another */
29
30@ The \CEE/ file \.{gb\_basic.c} has the following overall shape:
31
32@p
33#include "gb_graph.h" /* we use the {\sc GB\_\,GRAPH} data structures */
34@h@#
35@<Private variables@>@;
36@<Basic subroutines@>@;
37@<Applications of basic subroutines@>@;
38
39@ Several of the programs below allocate arrays that will be freed again
40before the routine is finished.
41
42@<Private variables@>=
43static Area working_storage;
44
45@ If a graph-generating subroutine encounters a problem, it returns |NULL|
46(that is, \.{NULL}), after putting a code number into the external variable
47|panic_code|. This code number identifies the type of failure.
48Otherwise the routine returns a pointer to the newly created graph, which
49will be represented with the data structures explained in {\sc GB\_\,GRAPH}.
50(The external variable |panic_code| is itself defined in {\sc GB\_\,GRAPH}.)
51
52@d panic(c)
53  {@+panic_code=c;
54    gb_free(working_storage);
55    gb_trouble_code=0;
56    return NULL;
57  }
58
59@ The names of vertices are sometimes formed from the names of other
60vertices, or from potentially long sequences of numbers. We assemble
61them in the |buffer| array, which is sufficiently long that the
62vast majority of applications will be unconstrained by size limitations.
63The programs assume that |BUF_SIZE| is rather large, but in cases of
64doubt they ensure that |BUF_SIZE| will never be exceeded.
65
66@d BUF_SIZE 4096
67
68@<Private v...@>=
69static char buffer[BUF_SIZE];
70
71@*Grids and game boards. The subroutine call
72|board(n1,n2,n3,n4,piece,wrap,directed)|
73constructs a graph based on the moves of generalized chesspieces on a
74generalized rectangular board. Each vertex of the graph corresponds to a
75position on the board. Each arc of the graph corresponds to a move from
76one position to another.
77
78The first parameters, |n1| through |n4|, specify the size of the board.
79If, for example, a two-dimensional board with $n_1$ rows and $n_2$ columns
80is desired, you set $|n1|=n_1$, $|n2|=n_2$, and $|n3|=0$; the resulting
81graph will have $n_1n_2$ vertices. If you want a three-dimensional
82board with $n_3$ layers, set $|n3|=n_3$ and $n_4=0$. If you want
83a 4-{\mc D} board, put the number of 4th coordinates in~|n4|.
84If you want a $d$-dimensional board with $2^d$ positions, set |n1=2|
85and |n2=-d|.
86
87In general, the |board| subroutine determines the dimensions by scanning the
88sequence |(n1,n2,n3,n4,0)=@t$(n_1,n_2,n_3,n_4,0)$@>| from left to right
89until coming to the first nonpositive parameter $n_{k+1}$. If $k=0$
90(i.e., if |n1<=0|), the default size $8\times8$ will be used; this is
91an ordinary chessboard with 8~rows and 8~columns. Otherwise if $n_{k+1}=0$,
92the board will have $k$~dimensions $n_1$, \dots,~$n_k$. Otherwise
93we must have $n_{k+1}<0$; in this case, the board will have $d=\vert n_{k+1}
94\vert$ dimensions, chosen as the first $d$ elements of the infinite
95periodic sequence $(n_1,\ldots,n_k,n_1,\ldots,n_k,n_1,\ldots\,)$.
96For example, the specification |(n1,n2,n3,n4)=(2,3,5,-7)| is about as
97tricky as you can get. It produces a seven-dimensional board with
98dimensions $(n_1,\ldots,n_7)=(2,3,5,2,3,5,2)$, hence a graph with
99$2\cdot3\cdot5\cdot2\cdot3\cdot5\cdot2=1800$ vertices.
100
101The |piece| parameter specifies the legal moves of a generalized chesspiece.
102If |piece>0|, a move from position~|u| to position~|v| is considered legal
103if and only if the Euclidean distance between points |u| and~|v| is
104equal to $\sqrt{\vphantom1\smash{|piece|}}$.
105For example, if |piece=1| and if we have a
106two-dimensional board, the legal moves from $(x,y)$ are to $(x,y\pm1)$ and
107$(x\pm1,y)$; these are the moves of a so-called wazir, the only moves that
108a king and a rook can both make. If |piece=2|, the legal moves from $(x,y)$
109are to $(x\pm1,y\pm1)$; these are the four moves that a king and a bishop
110can both make. (A piece that can make only these moves was called a ``fers''
111in ancient Muslim chess.) If |piece=5|, the legal moves are those of a
112knight, from $(x,y)$ to $(x\pm1,y\pm2)$ or to $(x\pm2,y\pm1)$. If |piece=3|,
113there are no legal moves on a two-dimensional board; but moves from
114$(x,y,z)$ to $(x\pm1,y\pm1,z\pm1)$ would be legal in three dimensions.
115If |piece=0|, it is changed to the default value |piece=1|.
116
117If the value of |piece| is negative, arbitrary multiples of the basic moves
118for $\vert|piece|\vert$ are permitted. For example, |piece=-1| defines the
119moves of a rook, from $(x,y)$ to $(x\pm a,y)$ or to $(x,y\pm a)$ for all
120$a>0$; |piece=-2| defines the moves of a bishop, from $(x,y)$ to
121$(x\pm a,y\pm a)$. The literature of ``fairy chess'' assigns standard names
122to the following |piece| values: $\rm wazir=1$, $\rm fers=2$, $\rm dabbaba=4$,
123$\rm knight=5$, $\rm alfil=8$, $\rm camel=10$, $\rm zebra=13$, $\rm giraffe
124=17$, $\rm fiveleaper=25$, $\hbox{root-50-leaper}=50$, etc.; $\rm rook=-1$,
125$\rm bishop=-2$, $\rm unicorn=-3$, $\rm dabbabarider=-4$, $\rm nightrider=-5$,
126$\rm alfilrider=-8$, $\rm camelrider=-10$, etc.
127
128To generate a board with the moves of a king, you can use the |gunion|
129subroutine below to take the union of boards with |piece=1| and
130|piece=2|. Similarly, you can get queen moves by taking the union of
131boards with |piece=-1| and |piece=-2|.
132
133If |piece>0|, all arcs of the graph will have length~1. If |piece<0|, the
134length of each arc will be the number of multiples of a basic move that
135produced the arc.
136
137@ If the |wrap| parameter is nonzero, it specifies a subset of coordinates
138in which values are computed modulo the corresponding size.
139For example, the coordinates $(x,y)$ for vertices on a two-dimensional
140board are restricted to the range $0\le x<n_1$, $0\le y<n_2$; therefore
141when |wrap=0|, a move from $(x,y)$ to $(x+\delta_1,y+\delta_2)$ is
142legal only if $0\le x+\delta_1<n_1$ and $0\le y+\delta_2<n_2$.
143But when |wrap=1|, the $x$~coordinates are allowed to ``wrap around'';
144the move would then be made to $((x+\delta_1)\bmod n_1,y+\delta_2)$,
145provided that $0\le y+\delta_2<n_2$. Setting |wrap=1| effectively
146makes the board into a cylinder instead of a rectangle. Similarly, the
147$y$~coordinates are allowed to wrap around when |wrap=2|. Both $x$ and
148$y$ coordinates are treated modulo their corresponding sizes when
149|wrap=3|; the board is then effectively a torus.  In general,
150coordinates $k_1$, $k_2$, \dots~ will wrap around when
151$|wrap|=2^{k_1-1}+2^{k_2-1}+\cdots\,$. Setting |wrap=-1| causes all
152coordinates to be computed modulo their size.
153
154The graph constructed by |board| will be undirected unless |directed!=0|.
155Directed |board| graphs will be acyclic when |wrap=0|, but they may
156have cycles when |wrap!=0|. Precise rules defining the directed arcs
157are given below.
158
159Several important special cases are worth noting. To get the complete graph
160on |n| vertices, you can say |board(n,0,0,0,-1,0,0)|. To get the
161transitive tournament on |n| vertices, i.e., the directed graph
162with arcs from |u| to |v| when |u<v|, you can say |board(n,0,0,0,-1,0,1)|.
163To get the empty graph on |n| vertices, you can say |board(n,0,0,0,2,0,0)|.
164To get a circuit (undirected) or a cycle (directed) of length~|n|,
165you can say |board(n,0,0,0,1,1,0)| and |board(n,0,0,0,1,1,1)|,
166respectively.
167
168@(gb_basic.h@>=
169#define complete(n) @[board((long)(n),0L,0L,0L,-1L,0L,0L)@]
170#define transitive(n) @[board((long)(n),0L,0L,0L,-1L,0L,1L)@]
171#define empty(n) @[board((long)(n),0L,0L,0L,2L,0L,0L)@]
172#define circuit(n) @[board((long)(n),0L,0L,0L,1L,1L,0L)@]
173#define cycle(n) @[board((long)(n),0L,0L,0L,1L,1L,1L)@]
174
175@ @<Basic subroutines@>=
176Graph *board(n1,n2,n3,n4,piece,wrap,directed)
177  long n1,n2,n3,n4; /* size of board desired */
178  long piece; /* type of moves desired */
179  long wrap; /* mask for coordinate positions that wrap around */
180  long directed; /* should the graph be directed? */
181{@+@<Vanilla local variables@>@;@+@q{@>@t}\6{@>@q}@>
182  long n; /* total number of vertices */
183  long p; /* $\vert|piece|\vert$ */
184  long l; /* length of current arc */
185  @<Normalize the board-size parameters@>;
186  @<Set up a graph with |n| vertices@>;
187  @<Insert arcs or edges for all legal moves@>;
188  if (gb_trouble_code) {
189    gb_recycle(new_graph);
190    panic(alloc_fault); /* alas, we ran out of memory somewhere back there */
191  }
192  return new_graph;
193}
194
195@ Most of the subroutines in {\sc GB\_\,BASIC} use the following local
196variables.
197
198@<Vanilla local variables@>=
199Graph *new_graph; /* the graph being constructed */
200register long i,j,k; /* all-purpose indices */
201register long d; /* the number of dimensions */
202register Vertex *v; /* the current vertex of interest */
203register long s; /* accumulator */
204
205@ Several arrays will facilitate the calculations that |board| needs to make.
206The number of distinct values in coordinate position~$k$ will be |nn[k]|;
207this coordinate position will wrap around if and only if |wr[k]!=0|.
208The current moves under consideration will be from $(x_1,\ldots,x_d)$
209to $(x_1+\delta_1,\ldots, x_d+\delta_d)$, where $\delta_k$ is stored
210in |del[k]|. An auxiliary array |sig| holds the sums
211$\sigma_k=\delta_1^2+\cdots+\delta_{k-1}^2$.  Additional arrays |xx|
212and |yy| hold coordinates of vertices before and after a move is made.
213
214Some of these arrays are also used for other purposes by other programs
215besides |board|; we will meet those programs later.
216
217We limit the number of dimensions to 91 or less. This is hardly a limitation,
218since the number of vertices would be astronomical even if the dimensionality
219were only half this big. But some of our later programs will be able
220to make good use of 40 or 50 dimensions and perhaps more; the number 91
221is an upper limit imposed by the number of standard printable characters
222(see the convention for vertex names in the |perms| routine).
223
224@d MAX_D 91
225
226@<Private...@>=
227static long nn[MAX_D+1]; /* component sizes */
228static long wr[MAX_D+1]; /* does this component wrap around? */
229static long del[MAX_D+1]; /* displacements for the current move */
230static long sig[MAX_D+2]; /* partial sums of squares of displacements */
231static long xx[MAX_D+1], yy[MAX_D+1]; /* coordinate values */
232
233@ @<Normalize the board-size parameters@>=
234if (piece==0) piece=1;
235if (n1<=0) {@+n1=n2=8;@+n3=0;@+}
236nn[1]=n1;
237if (n2<=0) {@+k=2;@+d=-n2;@+n3=n4=0;@+}
238else {
239  nn[2]=n2;
240  if (n3<=0) {@+k=3;@+d=-n3;@+n4=0;@+}
241  else {
242    nn[3]=n3;
243    if (n4<=0) {@+k=4;@+d=-n4;@+}
244    else {@+nn[4]=n4;@+d=4;@+goto done;@+}
245  }
246}
247if (d==0) {@+d=k-1;@+goto done;@+}
248@<Compute component sizes periodically for |d| dimensions@>;
249done: /* now |nn[1]| through |nn[d]| are set up */
250
251@ At this point, |nn[1]| through |nn[k-1]| are the component sizes
252that should be replicated periodically. In unusual cases, the number
253of dimensions might not be as large as the number of specifications.
254
255@<Compute component sizes periodically...@>=
256if (d>MAX_D) panic(bad_specs); /* too many dimensions */
257for (j=1; k<=d; j++,k++) nn[k]=nn[j];
258
259@ We want to make the subroutine idiot-proof, so we use floating-point
260arithmetic to make sure that boards with more than a billion cells have
261not been specified.
262
263@d MAX_NNN 1000000000.0
264
265@<Set up a graph with |n| vertices@>=
266{@+float nnn; /* approximate size */
267  for (n=1,nnn=1.0,j=1; j<=d; j++) {
268    nnn *= (float)nn[j];
269    if (nnn>MAX_NNN) panic(very_bad_specs); /* way too big */
270    n *= nn[j]; /* this multiplication cannot cause integer overflow */
271  }
272  new_graph=gb_new_graph(n);
273  if (new_graph==NULL)
274    panic(no_room); /* out of memory before we're even started */
275  sprintf(new_graph->id,"board(%ld,%ld,%ld,%ld,%ld,%ld,%d)",
276    n1,n2,n3,n4,piece,wrap,directed?1:0);
277  strcpy(new_graph->util_types,"ZZZIIIZZZZZZZZ");
278  @<Give names to the vertices@>;
279}
280
281@ The symbolic name of a board position like $(3,1)$ will be the string
282`\.{3.1}'. The first three coordinates are also stored as integers, in
283utility fields |x.I|, |y.I|, and |z.I|, because immediate access to
284those values will be helpful in certain applications. (The coordinates can,
285of course, always be recovered in a slower fashion from the vertex name,
286via |sscanf|.)
287
288The process of assigning coordinate values and names is equivalent to
289adding unity in a mixed-radix number system. Vertex $(x_1,\ldots,x_d)$
290will be in position $x_1n_2\ldots n_d+\cdots+x_{d-1}n_d+x_d$ relative
291to the first vertex of the new graph; therefore it is also possible to
292deduce the coordinates of a vertex from its address.
293
294@<Give names...@>=
295{@+register char *q; /* string pointer */
296  nn[0]=xx[0]=xx[1]=xx[2]=xx[3]=0;
297  for (k=4;k<=d;k++) xx[k]=0;
298  for (v=new_graph->vertices;;v++) {
299    q=buffer;
300    for (k=1;k<=d;k++) {
301      sprintf(q,".%ld",xx[k]);
302      while (*q) q++;
303    }
304    v->name=gb_save_string(&buffer[1]); /* omit |buffer[0]|, which is |'.'| */
305    v->x.I=xx[1];@+v->y.I=xx[2];@+v->z.I=xx[3];
306    for (k=d;xx[k]+1==nn[k];k--) xx[k]=0;
307    if (k==0) break; /* a ``carry'' has occurred all the way to the left */
308    xx[k]++; /* increase coordinate |k| */
309  }
310}
311
312@ Now we come to a slightly tricky part of the routine: the move generator.
313Let $p=\vert|piece|\vert$. The outer loop of this procedure runs through all
314solutions of the equation $\delta_1^2+\cdots+\delta_d^2=p$, where the
315$\delta$'s are nonnegative integers. Within that loop, we attach signs
316to the $\delta$'s, but we always leave $\delta_k$ positive if $\delta_1=
317\cdots=\delta_{k-1}=0$. For every such vector~$\delta$, we generate moves
318from |v| to $v+\delta$ for every vertex |v|. When |directed=0|,
319we use |gb_new_edge| instead of |gb_new_arc|, so that the reverse arc
320from $v+\delta$ to~|v| is also generated.
321
322@<Insert arcs or edges for all legal moves@>=
323@<Initialize the |wr|, |sig|, and |del| tables@>;
324p=piece;
325if (p<0) p=-p;
326while (1) {
327  @<Advance to the next nonnegative |del| vector, or |break| if done@>;
328  while (1) {
329    @<Generate moves for the current |del| vector@>;
330    @<Advance to the next signed |del| vector, or restore |del|
331          to nonnegative values and |break|@>;
332  }
333}
334
335@ The \CEE/ language does not define |>>| unambiguously. If |w| is negative,
336the assignment `|w>>=1|' here should keep |w| negative.
337(However, this technicality doesn't matter except in highly unusual cases
338when there are more than 32 dimensions.)
339@^system dependencies@>
340
341@<Initialize the |wr|, |sig|, and |del| tables@>=
342{@+register long w=wrap;
343  for (k=1;k<=d;k++,w>>=1) {
344    wr[k]=w&1;
345    del[k]=sig[k]=0;
346  }
347  sig[0]=del[0]=sig[d+1]=0;
348}
349
350@ The |sig| array makes it easy to backtrack through all partitions
351of |p| into an ordered sum of squares.
352
353@<Advance to the next nonnegative |del|...@>=
354for (k=d;sig[k]+(del[k]+1)*(del[k]+1)>p;k--) del[k]=0;
355if (k==0) break;
356del[k]++;
357sig[k+1]=sig[k]+del[k]*del[k];
358for (k++;k<=d;k++) sig[k+1]=sig[k];
359if (sig[d+1]<p) continue;
360
361@ @<Advance to the next signed |del| vector, or restore |del|
362          to nonnegative values and |break|@>=
363for (k=d;del[k]<=0;k--) del[k]=-del[k];
364if (sig[k]==0) break; /* all but |del[k]| were negative or zero */
365del[k]=-del[k]; /* some entry preceding |del[k]| is positive */
366
367@ We use the mixed-radix addition technique again when generating moves.
368
369@<Generate moves for the current |del| vector@>=
370for (k=1;k<=d;k++) xx[k]=0;
371for (v=new_graph->vertices;;v++) {
372  @<Generate moves from |v| corresponding to |del|@>;
373  for (k=d;xx[k]+1==nn[k];k--) xx[k]=0;
374  if (k==0) break; /* a ``carry'' has occurred all the way to the left */
375  xx[k]++; /* increase coordinate |k| */
376}
377
378@ The legal moves when |piece| is negative are derived as follows, in
379the presence of possible wraparound: Starting at $(x_1,\ldots,x_d)$, we
380move to $(x_1+\delta_1,\ldots,x_d+\delta_d)$, $(x_1+2\delta_1,\ldots,
381x_d+2\delta_d)$,~\dots, until either coming to a position with a nonwrapped
382coordinate out of range or coming back to the original point.
383
384A subtle technicality should be noted: When coordinates are wrapped and
385|piece>0|, self-loops are possible---for example, in |board(1,0,0,0,1,1,1)|.
386But self-loops never arise when |piece<0|.
387
388@<Generate moves from |v|...@>=
389for (k=1;k<=d;k++) yy[k]=xx[k]+del[k];
390for (l=1;;l++) {
391  @<Correct for wraparound, or |goto no_more| if off the board@>;
392  if (piece<0) @<Go to |no_more| if |yy=xx|@>;
393  @<Record a legal move from |xx| to |yy|@>;
394  if (piece>0) goto no_more;
395  for (k=1;k<=d;k++) yy[k]+=del[k];
396}
397no_more:@;
398
399@ @<Go to |no_more|...@>=
400{
401  for (k=1;k<=d;k++) if (yy[k]!=xx[k]) goto unequal;
402  goto no_more;
403  unequal:;
404}
405
406@ @<Correct for wraparound, or |goto no_more| if off the board@>=
407for (k=1;k<=d;k++) {
408  if (yy[k]<0) {
409    if (!wr[k]) goto no_more;
410    do yy[k]+=nn[k];@+ while (yy[k]<0);
411  }@+else if (yy[k]>=nn[k]) {
412    if (!wr[k]) goto no_more;
413    do yy[k]-=nn[k];@+ while (yy[k]>=nn[k]);
414  }
415}
416
417@ @<Record a legal move from |xx| to |yy|@>=
418for (k=2,j=yy[1];k<=d;k++) j=nn[k]*j+yy[k];
419if (directed) gb_new_arc(v,new_graph->vertices+j,l);
420else gb_new_edge(v,new_graph->vertices+j,l);
421
422@* Generalized triangular boards. The subroutine call
423|simplex(n,n0,n1,n2,n3,n4,directed)| creates a graph based on
424generalized triangular or tetrahedral configurations. Such graphs are
425similar in spirit to the game boards created by |board|, but they
426pertain to nonrectangular grids like those in Chinese checkers. As
427with |board| in the case |piece=1|, the vertices represent board positions
428and the arcs run from board positions to their nearest neighbors. Each arc has
429length~1.{\tolerance=1000\par}
430
431More formally, the vertices can be defined as sequences of nonnegative
432integers $(x_0,x_1,\ldots,x_d)$ whose sum is~|n|, where two sequences
433are considered adjacent if and only if they differ by $\pm1$ in exactly
434two components---equivalently, if the Euclidean distance between them
435is~$\sqrt2$. When $d=2$, for example, the vertices can be visualized
436as a triangular array
437$$\vcenter{\halign{&\hbox to 2em{\hss$#$\hss}\cr
438&&&(0,0,3)\cr
439&&(0,1,2)&&(1,0,2)\cr
440&(0,2,1)&&(1,1,1)&&(2,0,1)\cr
441(0,3,0)&&(1,2,0)&&(2,1,0)&&(3,0,0)\cr}}$$
442containing $(n+1)(n+2)/2$ elements, illustrated here when $n=3$; each vertex of
443the array has up to 6 neighbors. When $d=3$ the vertices form a tetrahedral
444array, a stack of triangular layers, and they can have as many as 12
445neighbors. In general, a vertex in a $d$-simplicial array will have up to
446$d(d+1)$ neighbors.
447
448If the |directed| parameter is nonzero, arcs run only from vertices to
449neighbors that are lexicographically greater---for example, downward
450or to the right in the triangular array shown. The directed graph is
451therefore acyclic, and a vertex of a $d$-simplicial array has
452out-degree at most $d(d+1)/2$.
453
454@ The first parameter, |n|, specifies the sum of the coordinates
455$(x_0,x_1,\ldots,x_d)$. The following parameters |n0| through |n4|
456specify upper bounds on those coordinates, and they also specify the
457dimensionality~|d|.
458
459If, for example, |n0|, |n1|, and |n2| are positive while |n3=0|, the
460value of~|d| will be~2 and the coordinates will be constrained to
461satisfy $0\le x_0\le|n0|$, $0\le x_1\le|n1|$, $0\le x_2\le|n2|$. These
462upper bounds essentially lop off the corners of the triangular array.
463We obtain a hexagonal board with $6m$ boundary cells by asking for
464|simplex(3m,2m,2m,2m,0,0,0)|. We obtain the diamond-shaped board used
465in the game of Hex [Martin Gardner, {\sl The Scientific American
466@^Gardner, Martin@>
467Book of Mathematical Puzzles {\char`\&} Diversions\/} (Simon {\char`\&}
468Schuster, 1959), Chapter~8] by calling |simplex(20,10,20,10,0,0,0)|.
469
470In general, |simplex| determines |d| and upper bounds $(n_0,n_1,\ldots,n_d)$
471in the following way: Let the first nonpositive entry of the sequence
472|(n0,n1,n2,n3,n4,0)|$\null=(n_0,n_1,n_2,n_3,n_4,0)$ be~$n_k$. If $k>0$
473and $n_k=0$, the value of~$d$ will be $k-1$ and the coordinates will be
474bounded by the given numbers $(n_0,\ldots,n_d)$. If $k>0$ and $n_k<0$,
475the value of~$d$ will be $\vert n_k\vert$ and the coordinates will be
476bounded by the first $d+1$ elements of the infinite periodic sequence
477$(n_0,\ldots,n_{k-1},n_0,\ldots,n_{k-1},n_0,\ldots\,)$. If $k=0$ and
478$n_0<0$, the value of~$d$ will be $\vert n_0\vert$ and the coordinates
479will be unbounded; equivalently, we may set $n_0=\cdots=n_d=n$. In
480this case the number of vertices will be $n+d\choose d$. Finally,
481if $k=0$ and $n_0=0$, we have the default case of a triangular array
482with $3n$ boundary cells, exactly as if $n_0=-2$.
483
484For example, the specification |n0=3|, |n1=-5| will produce all vertices
485$(x_0,x_1,\ldots,x_5)$ such that $x_0+x_1+\cdots+x_5=n$ and $0\le x_j\le3$.
486The specification |n0=1|, |n1=-d| will essentially produce all $n$-element
487subsets of the $(d+1)$-element set $\{0,1,\ldots,d\}$, because we can
488regard an element~$k$ as being present in the set if $x_k=1$ and absent
489if $x_k=0$. In that case two subsets are adjacent if and only if
490they have exactly $n-1$ elements in common.
491
492@ @<Basic subroutines@>=
493Graph *simplex(n,n0,n1,n2,n3,n4,directed)
494  unsigned long n; /* the constant sum of all coordinates */
495  long n0,n1,n2,n3,n4; /* constraints on coordinates */
496  long directed; /* should the graph be directed? */
497{@+@<Vanilla local variables@>@;@#
498  @<Normalize the simplex parameters@>;
499  @<Create a graph with one vertex for each point@>;
500  @<Name the points and create the arcs or edges@>;
501  if (gb_trouble_code) {
502    gb_recycle(new_graph);
503    panic(alloc_fault); /* darn, we ran out of memory somewhere back there */
504  }
505  return new_graph;
506}
507
508@ @<Normalize the simplex parameters@>=
509if (n0==0) n0=-2;
510if (n0<0) {@+k=2;@+nn[0]=n;@+d=-n0;@+n1=n2=n3=n4=0;@+}
511else {
512  if (n0>n) n0=n;
513  nn[0]=n0;
514  if (n1<=0) {@+k=2;@+d=-n1;@+n2=n3=n4=0;@+}
515  else {
516    if (n1>n) n1=n;
517    nn[1]=n1;
518    if (n2<=0) {@+k=3;@+d=-n2;@+n3=n4=0;@+}
519    else {
520      if (n2>n) n2=n;
521      nn[2]=n2;
522      if (n3<=0) {@+k=4;@+d=-n3;@+n4=0;@+}
523      else {
524        if (n3>n) n3=n;
525        nn[3]=n3;
526        if (n4<=0) {@+k=5;@+d=-n4;@+}
527        else {@+if (n4>n) n4=n;
528          nn[4]=n4;@+d=4;@+goto done;@+}
529      }
530    }
531  }
532}
533if (d==0) {@+d=k-2;@+goto done;@+}
534nn[k-1]=nn[0];
535@<Compute component sizes periodically...@>;
536done: /* now |nn[0]| through |nn[d]| are set up */
537
538@ @<Create a graph with one vertex for each point@>=
539@<Determine the number of feasible $(x_0,\ldots,x_d)$,
540  and allocate the graph@>;
541sprintf(new_graph->id,"simplex(%lu,%ld,%ld,%ld,%ld,%ld,%d)",
542    n,n0,n1,n2,n3,n4,directed?1:0);
543strcpy(new_graph->util_types,"VVZIIIZZZZZZZZ"); /* hash table will be used */
544
545@ We determine the number of vertices by determining the coefficient of~$z^n$
546in the power series
547$$(1+z+\cdots+z^{n_0})(1+z+\cdots+z^{n_1})\ldots(1+z+\cdots+z^{n_d}).$$
548
549@<Determine the number of feasible $(x_0,\ldots,x_d)$...@>=
550{@+long nverts; /* the number of vertices */
551  register long *coef=gb_typed_alloc(n+1,long,working_storage);
552  if (gb_trouble_code) panic(no_room+1); /* can't allocate |coef| array */
553  for (k=0;k<=nn[0];k++) coef[k]=1;
554   /* now |coef| represents the coefficients of $1+z+\cdots+z^{n_0}$ */
555  for (j=1;j<=d;j++)
556    @<Multiply the power series coefficients by $1+z+\cdots+z^{n_j}$@>;
557  nverts=coef[n];
558  gb_free(working_storage); /* recycle the |coef| array */
559  new_graph=gb_new_graph(nverts);
560  if (new_graph==NULL)
561    panic(no_room); /* out of memory before we're even started */
562}
563
564@ There's a neat way to multiply by $1+z+\cdots+z^{n_j}$: We multiply
565first by $1-z^{n_j+1}$, then sum the coefficients.
566
567We want to detect impossibly large specifications without risking
568integer overflow. It is easy to do this because multiplication is being
569done via addition.
570
571@<Multiply the power series coefficients by $1+z+\cdots+z^{n_j}$@>=
572{
573  for (k=n,i=n-nn[j]-1;i>=0;k--,i--) coef[k]-=coef[i];
574  s=1;
575  for (k=1;k<=n;k++) {
576    s+=coef[k];
577    if (s>1000000000) panic(very_bad_specs); /* way too big */
578    coef[k]=s;
579  }
580}
581
582@ As we generate the vertices, it proves convenient to precompute an
583array containing the numbers $y_j=n_j+\cdots+n_d$, which represent the
584largest possible sums $x_j+\cdots+x_d$. We also want to maintain
585the numbers $\sigma_j=n-(x_0+\cdots+x_{j-1})=x_j+\cdots+x_d$. The
586conditions
587$$0\le x_j\le n_j, \qquad \sigma_j-y_{j+1}\le x_j\le \sigma_j$$
588are necessary and sufficient, in the sense that we can find at least
589one way to complete a partial solution $(x_0,\ldots,x_k)$ to a full
590solution $(x_0,\ldots,x_d)$ if and only if the conditions hold for
591all $j\le k$.
592
593There is at least one solution if and only if $n\le y_0$.
594
595We enter the name string into a hash table, using the |hash_in|
596routine of {\sc GB\_\,GRAPH}, because there is no simple way to compute the
597location of a vertex from its coordinates.
598
599@<Name the points and create the arcs or edges@>=
600v=new_graph->vertices;
601yy[d+1]=0;@+sig[0]=n;
602for (k=d;k>=0;k--) yy[k]=yy[k+1]+nn[k];
603if (yy[0]>=n) {
604  k=0;@+xx[0]=(yy[1]>=n? 0: n-yy[1]);
605  while (1) {
606    @<Complete the partial solution $(x_0,\ldots,x_k)$@>;
607    @<Assign a symbolic name for $(x_0,\ldots,x_d)$ to vertex~|v|@>;
608    hash_in(v); /* enter |v->name| into the hash table
609                  (via utility fields |u,v|) */
610    @<Create arcs or edges from previous points to~|v|@>;
611    v++;
612    @<Advance to the next partial solution $(x_0,\ldots,x_k)$, where |k| is
613       as large as possible; |goto last| if there are no more solutions@>;
614  }
615}
616last:@+if (v!=new_graph->vertices+new_graph->n)
617  panic(impossible); /* can't happen */
618
619@ @<Complete the partial solution $(x_0,\ldots,x_k)$@>=
620for (s=sig[k]-xx[k],k++;k<=d;s-=xx[k],k++) {
621  sig[k]=s;
622  if (s<=yy[k+1]) xx[k]=0;
623  else xx[k]=s-yy[k+1];
624}
625if (s!=0) panic(impossible+1) /* can't happen */
626
627@ Here we seek the largest $k$ such that $x_k$ can be increased without
628violating the necessary and sufficient conditions stated earlier.
629
630@<Advance to the next partial solution $(x_0,\ldots,x_k)$...@>=
631for (k=d-1;;k--) {
632  if (xx[k]<sig[k] && xx[k]<nn[k]) break;
633  if (k==0) goto last;
634}
635xx[k]++;
636
637@ As in the |board| routine, we represent the sequence of coordinates
638$(2,0,1)$ by the string `\.{2.0.1}'.
639The string won't exceed |BUF_SIZE|, because the ratio |BUF_SIZE/MAX_D| is
640plenty big.
641
642The first three coordinate values, $(x_0,x_1,x_2)$, are placed into
643utility fields |x|, |y|, and |z|, so that they can be accessed immediately
644if an application needs them.
645
646@<Assign a symbolic name for $(x_0,\ldots,x_d)$ to vertex~|v|@>=
647{@+register char *p=buffer; /* string pointer */
648  for (k=0;k<=d;k++) {
649    sprintf(p,".%ld",xx[k]);
650    while (*p) p++;
651  }
652  v->name=gb_save_string(&buffer[1]); /* omit |buffer[0]|, which is |'.'| */
653  v->x.I=xx[0];@+v->y.I=xx[1];@+v->z.I=xx[2];
654}
655
656@ Since we are generating the vertices in lexicographic order of their
657coordinates, it is easy to identify all adjacent vertices that
658precede the current setting of $(x_0,x_1,\ldots,x_d)$. We locate them
659via their symbolic names.
660
661@<Create arcs or edges from previous points to~|v|@>=
662for (j=0;j<d;j++)
663  if (xx[j]) {@+register Vertex *u; /* previous vertex adjacent to |v| */
664    xx[j]--;
665    for (k=j+1;k<=d;k++)
666      if (xx[k]<nn[k]) {@+register char *p=buffer; /* string pointer */
667        xx[k]++;
668        for (i=0;i<=d;i++) {
669          sprintf(p,".%ld",xx[i]);
670          while (*p) p++;
671        }
672        u=hash_out(&buffer[1]);
673        if (u==NULL) panic(impossible+2); /* can't happen */
674        if (directed) gb_new_arc(u,v,1L);
675        else gb_new_edge(u,v,1L);
676        xx[k]--;
677      }
678    xx[j]++;
679  }
680
681@* Subset graphs. The subroutine call {\advance\thinmuskip 0mu plus 2mu
682|subsets(n,n0,n1,n2,n3,n4,size_bits,directed)|}
683creates a graph having the same vertices as
684|simplex(n,n0,n1,n2,n3,n4,directed)| but with a quite different notion
685of adjacency. In this we interpret a solution $(x_0,x_1,\ldots,x_d)$ to
686the conditions $x_0+x_1+\cdots+x_d=n$ and $0\le x_j\le n_j$ not as a
687position on a game board but as an $n$-element submultiset of the multiset
688$\{n_0\cdot0,n_1\cdot 1,\ldots,n_d\cdot d\}$ that has $x_j$ elements
689equal to~$j$. (If each $n_j=1$, the multiset is a set; this is an
690important special case.) Two vertices are adjacent if and only if
691their intersection has a cardinality that matches one of the bits in
692|size_bits|, which is an unsigned integer. Each arc has length~1.
693
694For example, suppose $n=3$ and |(n0,n1,n2,n3)=(2,2,2,0)|. Then the vertices
695are the 3-element submultisets of $\{0,0,1,1,2,2\}$, namely
696$$\{0,0,1\},\quad \{0,0,2\},\quad \{0,1,1\},\quad \{0,1,2\},\quad
697\{0,2,2\},\quad \{1,1,2\},\quad \{1,2,2\},$$
698which are represented by the respective vectors
699$$(2,1,0),\quad (2,0,1),\quad (1,2,0),\quad (1,1,1),\quad
700(1,0,2),\quad (0,2,1),\quad (0,1,2).$$
701The intersection of multisets represented by $(x_0,x_1,\ldots,x_d)$ and
702$(y_0,y_1,\ldots,y_d)$ is $$\bigl(\min(x_0,y_0),\min(x_1,y_1),\ldots,
703\min(x_d,y_d)\bigr);$$ each element occurs as often as it occurs
704in both multisets being intersected. If now |size_bits=3|, the
705multisets will be considered adjacent whenever their
706intersection contains exactly 0 or~1 elements, because $3=2^0+2^1$.
707The vertices adjacent to $\{0,0,1\}$, for example, will be
708$\{0,2,2\}$, $\{1,1,2\}$,
709 and $\{1,2,2\}$. In this case, every pair of submultisets
710has a nonempty intersection, so the same graph would be obtained
711if |size_bits=2|.
712
713If |directed| is nonzero, the graph will have directed arcs, from |u|
714to~|v| only if $u\le v$. Notice that the graph will have self-loops if
715and only if the binary representation of |size_bits| contains the term
716$2^n$, in which case there will be a loop from every vertex to itself.
717(In an undirected graph, such loops are represented by two arcs.)
718
719We define a macro |disjoint_subsets(n,k)| for the case
720of $n\choose k$ vertices, adjacent if and only if they represent
721disjoint $k$-subsets of an $n$-set.
722One important special case is the Petersen graph, whose vertices
723@^Petersen, Julius Peter Christian, graph@>
724are the 2-element subsets of $\{0,1,2,3,4\}$, adjacent when they
725are disjoint. This graph is remarkable because it contains 10 vertices,
726each of degree~3, but it has no circuits of length less than~5.
727
728@(gb_basic.h@>=
729#define disjoint_subsets(n,k) @[subsets((long)(k),1L,(long)(1-(n)),0L,0L,0L,1L,0L)@]
730#define petersen() @[disjoint_subsets(5,2)@]
731
732@ @<Basic subroutines@>=
733Graph *subsets(n,n0,n1,n2,n3,n4,size_bits,directed)
734  unsigned long n; /* the number of elements in the multiset */
735  long n0,n1,n2,n3,n4; /* multiplicities of elements */
736  unsigned long size_bits; /* intersection sizes that trigger arcs */
737  long directed; /* should the graph be directed? */
738{@+@<Vanilla local variables@>@;@#
739  @<Normalize the simplex parameters@>;
740  @<Create a graph with one vertex for each subset@>;
741  @<Name the subsets and create the arcs or edges@>;
742  if (gb_trouble_code) {
743    gb_recycle(new_graph);
744    panic(alloc_fault); /* rats, we ran out of memory somewhere back there */
745  }
746  return new_graph;
747}
748
749@ @<Create a graph with one vertex for each subset@>=
750@<Determine the number of feasible $(x_0,\ldots,x_d)$,
751  and allocate the graph@>;
752sprintf(new_graph->id,"subsets(%lu,%ld,%ld,%ld,%ld,%ld,0x%lx,%d)",
753    n,n0,n1,n2,n3,n4,size_bits,directed?1:0);
754strcpy(new_graph->util_types,"ZZZIIIZZZZZZZZ");
755 /* hash table will not be used */
756
757@ We generate the vertices with exactly the logic used in |simplex|.
758
759@<Name the subsets and create the arcs or edges@>=
760v=new_graph->vertices;
761yy[d+1]=0;@+sig[0]=n;
762for (k=d;k>=0;k--) yy[k]=yy[k+1]+nn[k];
763if (yy[0]>=n) {
764  k=0;@+xx[0]=(yy[1]>=n? 0: n-yy[1]);
765  while (1) {
766    @<Complete the partial solution $(x_0,\ldots,x_k)$@>;
767    @<Assign a symbolic name for $(x_0,\ldots,x_d)$ to vertex~|v|@>;
768    @<Create arcs or edges from previous subsets to~|v|@>;
769    v++;
770    @<Advance to the next partial solution $(x_0,\ldots,x_k)$, where |k| is
771       as large as possible; |goto last| if there are no more solutions@>;
772  }
773}
774last:@+if (v!=new_graph->vertices+new_graph->n)
775  panic(impossible); /* can't happen */
776
777@ The only difference is that we generate the arcs or edges by brute
778force, examining each pair of vertices to see if they are adjacent or not.
779
780The code here is character-set dependent: It assumes that `\..' and null
781have a character code less than `\.0', as in ASCII. It also assumes
782that characters occupy exactly eight bits.
783@^system dependencies@>
784
785@d UL_BITS 8*sizeof(unsigned long) /* the number of bits in |size_bits| */
786
787@<Create arcs or edges from previous subsets to~|v|@>=
788{@+register Vertex *u;
789  for (u=new_graph->vertices;u<=v;u++) {@+register char *p=u->name;
790    long ss=0; /* the number of elements common to |u| and |v| */
791    for (j=0;j<=d;j++,p++) {
792      for (s=(*p++)-'0';*p>='0';p++) s=10*s+*p-'0'; /* |sscanf(p,"%ld",&s)| */
793@^character-set dependencies@>
794      if (xx[j]<s) ss+=xx[j];
795      else ss+=s;
796    }
797    if (ss<UL_BITS && (size_bits&(((unsigned long)1)<<ss))) {
798      if (directed) gb_new_arc(u,v,1L);
799      else gb_new_edge(u,v,1L);
800    }
801  }
802}
803
804@* Permutation graphs. The subroutine call
805|perms(n0,n1,n2,n3,n4,max_inv,directed)|
806creates a graph whose vertices represent the permutations of a
807multiset that have at most |max_inv| inversions. Two permutations are adjacent
808in the graph if one is obtained from the other by interchanging two
809adjacent elements. Each arc has length~1.
810
811For example, the multiset $\{0,0,1,2\}$ has the following twelve permutations:
812$$\vcenter{\halign{#&&\quad#\cr
8130012,&0021,&0102,&0120,&0201,&0210,\cr
8141002,&1020,&1200,&2001,&2010,&2100.\cr}}$$
815The first of these, 0012, has two neighbors, 0021 and 0102.
816
817The number of inversions is the number of pairs of elements $xy$ such
818that $x>y$ and $x$ precedes $y$ from left to right, counting
819multiplicity.  For example, 2010 has four inversions, corresponding to
820$xy\in\{20,21,20,10\}$. It is not difficult to verify that the number
821of inversions of a permutation is the distance in the graph from that
822permutation to the lexicographically first permutation.
823
824Parameters |n0| through |n4| specify the composition of the multiset,
825just as in the |subsets| routine.
826Roughly speaking, there are |n0| elements equal to~0, |n1| elements
827equal to~1, and so on. The multiset $\{0,0,1,2,3,3\}$, for example, would
828be represented by |(n0,n1,n2,n3,n4)=(2,1,1,2,0)|.
829
830Of course, we sometimes want to have multisets with more than five distinct
831elements; when there are $d+1$ distinct elements, the multiset should have
832$n_k$ elements equal to~$k$ and $n=n_0+n_1+\cdots+n_d$ elements in all.
833Larger values of $d$ can be specified by using |-d| as a parameter:
834If |n0=-d|, each multiplicity $n_k$ is taken to be~1; if |n0>0| and |n1=-d|,
835each multiplicity $n_k$ is taken to be equal to~|n0|;
836if |n0>0|, |n1>0|, and |n2=-d|, the multiplicities are alternately
837$(|n0|,|n1|,|n0|,|n1|,|n0|,\ldots\,)$; if  |n0>0|, |n1>0|, |n2>0|, and
838|n3=-d|, the multiplicities are the first~|d+1| elements of the
839periodic sequence $(|n0|,|n1|,|n2|,|n0|,|n1|,\ldots\,)$; and if all
840but |n4| are positive, while |n4=-d|, the multiplicities again
841are periodic.
842
843An example like |(n0,n1,n2,n3,n4)=(1,2,3,4,-8)| is about as tricky
844as you can get. It specifies the multiset $\{0,1,1,2,2,2,3,3,3,3,4,5,5,
8456,6,6,7,7,7,7,8\}$.
846
847If any of the multiplicity parameters is negative or zero, the
848remaining multiplicities are ignored. For example, if |n2<=0|, the
849subroutine does not look at |n3| or~|n4|.
850
851You probably don't want to try |perms(n0,0,0,0,0,max_inv,directed)|
852when |n0>0|, because a multiset with |n0| identical elements has only
853one permutation.
854
855The special case when you want all $n!$ permutations of an $n$-element set
856can be obtained by calling |all_perms(n,directed)|.
857
858@(gb_basic.h@>=
859#define all_perms(n,directed) @[perms((long)(1-(n)),0L,0L,0L,0L,0L,\
860   (long)(directed))@]
861
862@ If |max_inv=0|, all permutations will be considered, regardless of
863the number of inversions. In that case the total number of vertices in
864the graph will be the multinomial coefficient $${n\choose
865n_0,n_1,\ldots,n_d}\,,\qquad n=n_0+n_1+\cdots+n_d.$$ The maximum
866number of inversions in general is the number of inversions of the
867lexicographically last permutation, namely ${n\choose2}-{n_0\choose2}-
868{n_1\choose2}-\cdots-{n_d\choose2}=\sum_{0\le j<k\le d}n_jn_k$.
869
870\vskip1pt
871Notice that in the case $d=1$, we essentially obtain all combinations of
872|n0+n1| elements taken |n1| at a time. The positions of the 1's correspond
873to the elements of a subset or sample.
874
875If |directed| is nonzero, the graph will contain only directed arcs
876from permutations to neighboring permutations that have exactly one more
877inversion. In this case the graph corresponds to a partial ordering that is a
878lattice with interesting properties; see the article by Bennett and Birkhoff
879@^Bennett, Mary Katherine@>
880@^Birkhoff, Garrett@>
881in {\sl Algebra Universalis\/ \bf32} (1994), 115--144.
882
883@ The program for |perms| is very similar in structure to the program
884for |simplex| already considered.
885
886@<Basic subroutines@>=
887Graph *perms(n0,n1,n2,n3,n4,max_inv,directed)
888  long n0,n1,n2,n3,n4; /* composition of the multiset */
889  unsigned long max_inv; /* maximum number of inversions */
890  long directed; /* should the graph be directed? */
891{@+@<Vanilla local variables@>@;@+@q{@>@t}\6{@>@q}@>
892  register long n; /* total number of elements in multiset */
893  @<Normalize the permutation parameters@>;
894  @<Create a graph with one vertex for each permutation@>;
895  @<Name the permutations and create the arcs or edges@>;
896  if (gb_trouble_code) {
897    gb_recycle(new_graph);
898    panic(alloc_fault); /* shucks, we ran out of memory somewhere back there */
899  }
900  return new_graph;
901}
902
903@ @<Normalize the permutation parameters@>=
904if (n0==0) {@+n0=1;@+n1=0;@+} /* convert the empty set into $\{0\}$ */
905else if (n0<0) {@+n1=n0;@+n0=1;@+}
906n=BUF_SIZE; /* this allows us to borrow code from |simplex|, already written */
907@<Normalize the simplex parameters@>;
908@<Determine |n| and the maximum possible number of inversions@>;
909
910@ Here we want to set |max_inv| to the maximum possible number of
911inversions, if the given value of |max_inv| is zero or
912if it exceeds that maximum number.
913
914@<Determine |n| and the maximum possible number of inversions@>=
915{@+register long ss; /* max inversions known to be possible */
916  for (k=0,s=ss=0;k<=d;ss+=s*nn[k],s+=nn[k],k++)
917    if (nn[k]>=BUF_SIZE) panic(bad_specs);
918      /* too many elements in the multiset */
919  if (s>=BUF_SIZE) panic(bad_specs+1); /* too many elements in the multiset */
920  n=s;
921  if (max_inv==0 || max_inv>ss) max_inv=ss;
922}
923
924@ To determine the number of vertices, we sum the first |max_inv+1|
925coefficients of a power series in which the coefficient of~$z^j$
926is the number of permutations having $j$ inversions. It is known
927[{\sl Sorting and Searching}, exercise 5.1.2--16] that this power series
928is the ``$z$-multinomial coefficient''
929$${n\choose n_0,\ldots,n_d}_{\!z}={n!_z\over n_0!_z\ldots n_d!_z}\,,
930\qquad\hbox{where}\qquad m!_z=\prod_{k=1}^m{1-z^k\over 1-z}\,.$$
931
932@<Create a graph with one vertex for each permutation@>=
933{@+long nverts; /* the number of vertices */
934  register long *coef=gb_typed_alloc(max_inv+1,long,working_storage);
935  if (gb_trouble_code) panic(no_room+1); /* can't allocate |coef| array */
936  coef[0]=1;
937  for (j=1,s=nn[0];j<=d;s+=nn[j],j++)
938    @<Multiply the power series coefficients by
939       $\prod_{1\le k\le n_j}(1-z^{s+k})/(1-z^k)$@>;
940  for (k=1,nverts=1;k<=max_inv;k++) {
941    nverts+=coef[k];
942    if (nverts>1000000000) panic(very_bad_specs); /* way too big */
943  }
944  gb_free(working_storage); /* recycle the |coef| array */
945  new_graph=gb_new_graph(nverts);
946  if (new_graph==NULL)
947    panic(no_room); /* out of memory before we're even started */
948  sprintf(new_graph->id,"perms(%ld,%ld,%ld,%ld,%ld,%lu,%d)",
949    n0,n1,n2,n3,n4,max_inv,directed?1:0);
950  strcpy(new_graph->util_types,"VVZZZZZZZZZZZZ"); /* hash table will be used */
951}
952
953@ After multiplication by $(1-z^{k+s})/(1-z^k)$, the coefficients of the
954power series will be nonnegative, because they are the coefficients of
955a $z$-multinomial coefficient.
956
957@<Multiply the power series coefficients by
958       $\prod_{1\le k\le n_j}(1-z^{s+k})/(1-z^k)$@>=
959for (k=1;k<=nn[j];k++) {@+register long ii;
960  for (i=max_inv,ii=i-k-s;ii>=0;ii--,i--) coef[i]-=coef[ii];
961  for (i=k,ii=0;i<=max_inv;i++,ii++) {
962    coef[i]+=coef[ii];
963    if (coef[i]>1000000000) panic(very_bad_specs+1); /* way too big */
964  }
965}
966
967@ As we generate the permutations, we maintain a table $(y_1,\ldots,y_n)$,
968where $y_k$ is the number of
969inversions whose first element is the $k$th element of the multiset.
970For example, if the multiset is $\{0,0,1,2\}$ and the current permutation is
971$(2,0,1,0)$, the inversion table is $(y_1,y_2,y_3,y_4)=(0,0,1,3)$. Clearly
972$0\le y_k<k$, and $y_k\le y_{k-1}$ when the $k$th element of the multiset
973is the same as the $(k-1)$st element. These conditions are necessary
974and sufficient to define a valid inversion table. We will generate
975permutations in lexicographic order of their inversion tables.
976
977For convenience, we set up another array~|z|, which holds the
978initial inversion-free permutation.
979
980@<Name the permutations and create the arcs or edges@>=
981{@+register long *xtab,*ytab,*ztab; /* permutations and their inversions */
982  long m=0; /* current number of inversions */
983  @<Initialize |xtab|, |ytab|, and |ztab|@>;
984  v=new_graph->vertices;
985  while (1) {
986    @<Assign a symbolic name for $(x_1,\ldots,x_n)$ to vertex~|v|@>;
987    @<Create arcs or edges from previous permutations to~|v|@>;
988    v++;
989    @<Advance to the next perm; |goto last| if there are no more solutions@>;
990  }
991  last:@+if (v!=new_graph->vertices+new_graph->n)
992    panic(impossible); /* can't happen */
993  gb_free(working_storage);
994}
995
996@ @<Initialize |xtab|, |ytab|, and |ztab|@>=
997xtab=gb_typed_alloc(3*n+3,long,working_storage);
998if (gb_trouble_code) { /* can't allocate |xtab| */
999 gb_recycle(new_graph);@+panic(no_room+2);@+}
1000ytab=xtab+(n+1);
1001ztab=ytab+(n+1);
1002for (j=0,k=1,s=nn[0];;k++) {
1003  xtab[k]=ztab[k]=j; /* |ytab[k]=0| */
1004  if (k==s) {
1005    if (++j>d) break;
1006    else s+=nn[j];
1007  }
1008}
1009
1010@ Here is the heart of the permutation logic. We find the largest~$k$
1011such that $y_k$ can legitimately be increased by~1. When we encounter
1012a~$k$ for which $y_k$ cannot be increased, we set $y_k=0$ and adjust
1013the $x$'s accordingly. If no $y_k$ can be increased, we are done.
1014
1015@<Advance to the next perm...@>=
1016for (k=n;k;k--) {
1017  if (m<max_inv && ytab[k]<k-1)
1018    if (ytab[k]<ytab[k-1] || ztab[k]>ztab[k-1]) goto move;
1019  if (ytab[k]) {
1020    for (j=k-ytab[k];j<k;j++) xtab[j]=xtab[j+1];
1021    m-=ytab[k];
1022    ytab[k]=0;
1023    xtab[k]=ztab[k];
1024  }
1025}
1026goto last;
1027move: j=k-ytab[k]; /* the current location of the $k$th element, $z_k$ */
1028xtab[j]=xtab[j-1];@+xtab[j-1]=ztab[k];
1029ytab[k]++;@+m++;
1030
1031@ A permutation is encoded as a sequence of nonblank characters,
1032using an abbreviated copy of the |imap| code from {\sc GB\_\,IO} and omitting
1033the characters that need to be quoted within strings. If the
1034number of distinct elements in the multiset is at most~62, only digits
1035and letters will appear in the vertex name.
1036
1037@<Private variables@>=
1038static char *short_imap="0123456789ABCDEFGHIJKLMNOPQRSTUVWXYZ\
1039abcdefghijklmnopqrstuvwxyz_^~&@@,;.:?!%#$+-*/|<=>()[]{}`'";
1040
1041@ @<Assign a symbolic name for $(x_1,\ldots,x_n)...@>=
1042{@+register char *p; register long *q;
1043  for (p=&buffer[n-1],q=&xtab[n];q>xtab;p--,q--) *p=short_imap[*q];
1044  v->name=gb_save_string(buffer);
1045  hash_in(v); /* enter |v->name| into the hash table
1046                  (via utility fields |u,v|) */
1047}
1048
1049@ Since we are generating the vertices in lexicographic order of their
1050inversions, it is easy to identify all adjacent vertices that
1051precede the current setting of $(x_1,\ldots,x_n)$. We locate them
1052via their symbolic names.
1053
1054@<Create arcs or edges from previous permutations to~|v|@>=
1055for (j=1;j<n;j++)
1056  if (xtab[j]>xtab[j+1]) {@+register Vertex *u;
1057                         /* previous vertex adjacent to |v| */
1058    buffer[j-1]=short_imap[xtab[j+1]];@+buffer[j]=short_imap[xtab[j]];
1059    u=hash_out(buffer);
1060    if (u==NULL) panic(impossible+2); /* can't happen */
1061    if (directed) gb_new_arc(u,v,1L);
1062    else gb_new_edge(u,v,1L);
1063    buffer[j-1]=short_imap[xtab[j]];@+buffer[j]=short_imap[xtab[j+1]];
1064  }
1065
1066@* Partition graphs. The subroutine call
1067|parts(n,max_parts,max_size,directed)|
1068creates a graph whose vertices represent the different ways to partition
1069the integer~|n| into at most |max_parts| parts, where each part is at most
1070|max_size|. Two partitions are adjacent in the graph if
1071one can be obtained from the other by combining two parts.
1072Each arc has length~1.
1073
1074For example, the partitions of~5 are
1075$$5,\quad 4+1,\quad 3+2,\quad 3+1+1,\quad 2+2+1,
1076       \quad 2+1+1+1,\quad 1+1+1+1+1.$$
1077Here 5 is adjacent to $4+1$ and to $3+2$; $4+1$ is adjacent also to
1078$3+1+1$ and to $2+2+1$; $3+2$ is adjacent also to $3+1+1$ and to $2+2+1$; etc.
1079If |max_size| is 3, the partitions 5 and $4+1$ would not be included in
1080the graph. If |max_parts| is 3, the partitions $2+1+1+1$ and $1+1+1+1+1$
1081would not be included.
1082
1083If |max_parts| or |max_size| are zero, they are reset to be equal to~|n|
1084so that they make no restriction on the partitions.
1085
1086If |directed| is nonzero, the graph will contain only directed arcs from
1087partitions to their neighbors that have exactly one more part.
1088
1089The special case when we want to generate all $p(n)$ partitions of the
1090integer~$n$ can be obtained by calling |all_parts(n,directed)|.
1091
1092@(gb_basic.h@>=
1093#define all_parts(n,directed) @[parts((long)(n),0L,0L,(long)(directed))@]
1094
1095@ The program for |parts| is very similar in structure to the program
1096for |perms| already considered.
1097
1098@<Basic subroutines@>=
1099Graph *parts(n,max_parts,max_size,directed)
1100  unsigned long n; /* the number being partitioned */
1101  unsigned long max_parts; /* maximum number of parts */
1102  unsigned long max_size; /* maximum size of each part */
1103  long directed; /* should the graph be directed? */
1104{@+@<Vanilla local variables@>@;@#
1105  if (max_parts==0 || max_parts>n) max_parts=n;
1106  if (max_size==0 || max_size>n) max_size=n;
1107  if (max_parts>MAX_D) panic(bad_specs); /* too many parts allowed */
1108  @<Create a graph with one vertex for each partition@>;
1109  @<Name the partitions and create the arcs or edges@>;
1110  if (gb_trouble_code) {
1111    gb_recycle(new_graph);
1112    panic(alloc_fault);
1113     /* doggone it, we ran out of memory somewhere back there */
1114  }
1115  return new_graph;
1116}
1117
1118@ The number of vertices is the coefficient of $z^n$
1119in the $z$-binomial coefficient ${m+p\choose m}_z$, where $m=|max_parts|$
1120and $p=|max_size|$. This coefficient is calculated as in the |perms| routine.
1121
1122@<Create a graph with one vertex for each partition@>=
1123{@+long nverts; /* the number of vertices */
1124  register long *coef=gb_typed_alloc(n+1,long,working_storage);
1125  if (gb_trouble_code) panic(no_room+1); /* can't allocate |coef| array */
1126  coef[0]=1;
1127  for (k=1;k<=max_parts;k++) {
1128    for (j=n,i=n-k-max_size;i>=0;i--,j--) coef[j]-=coef[i];
1129    for (j=k,i=0;j<=n;i++,j++) {
1130      coef[j]+=coef[i];
1131      if (coef[j]>1000000000) panic(very_bad_specs); /* way too big */
1132    }
1133  }
1134  nverts=coef[n];
1135  gb_free(working_storage); /* recycle the |coef| array */
1136  new_graph=gb_new_graph(nverts);
1137  if (new_graph==NULL)
1138    panic(no_room); /* out of memory before we're even started */
1139  sprintf(new_graph->id,"parts(%lu,%lu,%lu,%d)",
1140    n,max_parts,max_size,directed?1:0);
1141  strcpy(new_graph->util_types,"VVZZZZZZZZZZZZ"); /* hash table will be used */
1142}
1143
1144@ As we generate the partitions, we maintain
1145the numbers $\sigma_j=n-(x_1+\cdots+x_{j-1})=x_j+x_{j+1}+\cdots\,$,
1146somewhat as we did in the |simplex| routine. We set $x_0=|max_size|$
1147and $y_j=|max_parts|+1-j$; then when values $(x_1,\ldots,x_{j-1})$ are
1148given, the conditions
1149$$\sigma_j/y_j\le x_j\le \sigma_j,\qquad x_j\le x_{j-1}$$
1150characterize the legal values of~$x_j$.
1151
1152@<Name the partitions and create the arcs or edges@>=
1153v=new_graph->vertices;
1154xx[0]=max_size;@+sig[1]=n;
1155for (k=max_parts,s=1;k>0;k--,s++) yy[k]=s;
1156if (max_size*max_parts>=n) {
1157  k=1;@+xx[1]=(n-1)/max_parts+1; /* $\lceil n/|max_parts|\rceil$ */
1158  while (1) {
1159    @<Complete the partial solution $(x_1,\ldots,x_k)$@>;
1160    @<Assign the name $x_1+\cdots+x_d$ to vertex~|v|@>;
1161    @<Create arcs or edges from |v| to previous partitions@>;
1162    v++;
1163    @<Advance to the next partial solution $(x_1,\ldots,x_k)$, where |k| is
1164       as large as possible; |goto last| if there are no more solutions@>;
1165  }
1166}
1167last:@+if (v!=new_graph->vertices+new_graph->n)
1168  panic(impossible); /* can't happen */
1169
1170@ @<Complete the partial solution $(x_1,\ldots,x_k)$@>=
1171for (s=sig[k]-xx[k],k++;s;k++) {
1172  sig[k]=s;
1173  xx[k]=(s-1)/yy[k]+1;
1174  s-=xx[k];
1175}
1176d=k-1; /* the smallest part is $x_d$ */
1177
1178@ Here we seek the largest $k$ such that $x_k$ can be increased without
1179violating the necessary and sufficient conditions stated earlier.
1180
1181@<Advance to the next partial solution $(x_1,\ldots,x_k)$...@>=
1182if (d==1) goto last;
1183for (k=d-1;;k--) {
1184  if (xx[k]<sig[k] && xx[k]<xx[k-1]) break;
1185  if (k==1) goto last;
1186}
1187xx[k]++;
1188
1189@ @<Assign the name $x_1+...@>=
1190{@+register char *p=buffer; /* string pointer */
1191  for (k=1;k<=d;k++) {
1192    sprintf(p,"+%ld",xx[k]);
1193    while (*p) p++;
1194  }
1195  v->name=gb_save_string(&buffer[1]); /* omit |buffer[0]|, which is |'+'| */
1196  hash_in(v); /* enter |v->name| into the hash table
1197                  (via utility fields |u,v|) */
1198}
1199
1200@ Since we are generating the partitions in lexicographic order of their
1201parts, it is reasonably easy to identify all adjacent vertices that
1202precede the current setting of $(x_1,\ldots,x_d)$, by splitting
1203$x_j$ into two parts when $x_j\ne x_{j+1}$. We locate previous partitions
1204via their symbolic names.
1205
1206@<Create arcs or edges from |v| to previous partitions@>=
1207if (d<max_parts) {
1208  xx[d+1]=0;
1209  for (j=1;j<=d;j++) {
1210    if (xx[j]!=xx[j+1]) {@+long a,b;
1211      for (b=xx[j]/2,a=xx[j]-b;b;a++,b--)
1212        @<Generate a subpartition $(n_1,\ldots,n_{d+1})$ by
1213              splitting $x_j$ into $a+b$, and make that subpartition
1214              adjacent to~|v|@>;
1215    }
1216    nn[j]=xx[j];
1217  }
1218}
1219
1220@ The values of $(x_1,\ldots,x_{j-1})$ have already been copied into
1221$(n_1,\ldots,n_{j-1})$. Our job is to copy the smaller parts
1222$(x_{j+1},\ldots,x_d)$ while
1223inserting $a$ and $b$ in their proper places, knowing that $a\ge b$.
1224
1225@<Generate a subpartition $(n_1,\ldots,n_{d+1})$...@>=
1226{@+register Vertex *u; /* previous vertex adjacent to |v| */
1227  register char *p=buffer;
1228  for (k=j+1;xx[k]>a;k++) nn[k-1]=xx[k];
1229  nn[k-1]=a;
1230  for (;xx[k]>b;k++) nn[k]=xx[k];
1231  nn[k]=b;
1232  for (;k<=d;k++) nn[k+1]=xx[k];
1233  for (k=1;k<=d+1;k++) {
1234    sprintf(p,"+%ld",nn[k]);
1235    while (*p) p++;
1236  }
1237  u=hash_out(&buffer[1]);
1238  if (u==NULL) panic(impossible+2); /* can't happen */
1239  if (directed) gb_new_arc(v,u,1L);
1240  else gb_new_edge(v,u,1L);
1241}
1242
1243@* Binary tree graphs. The subroutine call
1244|binary(n,max_height,directed)|
1245creates a graph whose vertices represent the binary trees with $n$ internal
1246nodes and with all leaves at a distance that is at most |max_height|
1247from the root. Two binary trees are adjacent in the graph if
1248one can be obtained from the other by a single application of the
1249associative law for binary operations, i.e., by replacing some subtree
1250of the form $(\alpha\cdot\beta)\cdot\gamma$ by the subtree $\alpha\cdot
1251(\beta\cdot\gamma)$. (This transformation on binary trees is often
1252called a ``rotation.'') If the |directed| parameter is nonzero, the
1253directed arcs go from a tree containing $(\alpha\cdot\beta)\cdot\gamma$
1254to a tree containing $\alpha\cdot(\beta\cdot\gamma)$ in its place; otherwise
1255the graph is undirected. Each arc has length~1.
1256
1257For example, the binary trees with three internal nodes form a circuit of
1258length~5. They are
1259$$\mathcode`.="2201 % \cdot
1260(a.b).(c.d),\quad a.(b.(c.d)),\quad a.((b.c).d),\quad (a.(b.c)).d,\quad
1261((a.b).c).d,$$
1262if we use infix notation and name the leaves $(a,b,c,d)$ from left to right.
1263Here each tree is related to its two neighbors by associativity. The
1264first and last trees are also related in the same way.
1265
1266If |max_height=0|, it is changed to |n|, which means there is no
1267restriction on the height of a leaf. In this case the graph will have
1268exactly ${2n+1\choose n}/ (2n+1)$ vertices; furthermore, each vertex
1269will have exactly $n-1$ neighbors, because a rotation will be possible
1270just above every internal node except the root.  The graph in this
1271case can also be interpreted geometrically: The vertices are in
1272one-to-one correspondence with the triangulations of a regular
1273$(n+2)$-gon; two triangulations are adjacent if and only if one is obtained
1274from the other by replacing the pair of adjacent triangles $ABC,DCB$
1275by the pair $ADC,BDA$.
1276
1277The partial ordering corresponding to the directed graph on
1278${2n+1\choose n}/(2n+1)$ vertices created by |all_trees(n,1)|
1279is a lattice with interesting properties. See Huang and Tamari,
1280@^Huang, Samuel Shung@>
1281@^Tamari, Dov@>
1282{\sl Journal of Combinatorial Theory\/ \bf A13} (1972), 7--13.
1283
1284@(gb_basic.h@>=
1285#define all_trees(n,directed) @[binary((long)(n),0L,(long)(directed))@]
1286
1287@ The program for |binary| is very similar in structure to the program
1288for |parts| already considered. But the details are more exciting.
1289
1290@<Basic subroutines@>=
1291Graph *binary(n,max_height,directed)
1292  unsigned long n; /* the number of internal nodes */
1293  unsigned long max_height; /* maximum height of a leaf */
1294  long directed; /* should the graph be directed? */
1295{@+@<Vanilla local variables@>@;@#
1296  if (2*n+2>BUF_SIZE) panic(bad_specs); /* |n| is too huge for us */
1297  if (max_height==0 || max_height>n) max_height=n;
1298  if (max_height>30) panic(very_bad_specs); /* more than a billion vertices */
1299  @<Create a graph with one vertex for each binary tree@>;
1300  @<Name the trees and create the arcs or edges@>;
1301  if (gb_trouble_code) {
1302    gb_recycle(new_graph);
1303    panic(alloc_fault); /* uff da, we ran out of memory somewhere back there */
1304  }
1305  return new_graph;
1306}
1307
1308@ The number of vertices is the coefficient of $z^n$
1309in the power series $G_h$, where $h=|max_height|$ and $G_h$ satisfies
1310the recurrence
1311$$G_0=1,\qquad G_{h+1}=1+z G_h^2.$$
1312
1313The coefficients of $G_5$ are $\le55308$, but the
1314coefficients of $G_6$ are much larger; they exceed one billion when
1315$28\le n\le49$, and they exceed one million when $17\le n\le 56$.
1316In order to avoid overflow during this calculation, we use a
1317special method when $h\ge6$ and $n\ge20$: In such cases, graphs
1318of reasonable size arise only if $n\ge 2^h-7$, and we look at the
1319coefficient of $z^{-(2^h-1-n)}$ in $R_h=G_h/z^{2^h-1}$, which is a
1320power series in $z^{-1}$ defined by the recurrence
1321$$R_0=1,\qquad R_{h+1}=R_h^2+z^{1-2^{h+1}}.$$
1322
1323@<Create a graph with one vertex for each binary tree@>=
1324{@+long nverts; /* the number of vertices */
1325  if (n>=20 && max_height>=6) @<Compute |nverts| using the $R$ series@>@;
1326  else {
1327    nn[0]=nn[1]=1;
1328    for (k=2;k<=n;k++) nn[k]=0;
1329    for (j=2;j<=max_height;j++)
1330      for (k=n-1;k;k--) {
1331        for (s=0,i=k;i>=0;i--) s+=nn[i]*nn[k-i]; /* overflow impossible */
1332        nn[k+1]=s;
1333      }
1334    nverts=nn[n];
1335  }
1336  new_graph=gb_new_graph(nverts);
1337  if (new_graph==NULL)
1338    panic(no_room); /* out of memory before we're even started */
1339  sprintf(new_graph->id,"binary(%lu,%lu,%d)",
1340    n,max_height,directed?1:0);
1341  strcpy(new_graph->util_types,"VVZZZZZZZZZZZZ"); /* hash table will be used */
1342}
1343
1344@ The smallest nontrivial graph that is unilaterally disallowed by
1345this procedure on the grounds of size limitations occurs when |max_height=6|
1346and |n=20|; it has 14,162,220 vertices.
1347
1348@<Compute |nverts| using the $R$ series@>=
1349{@+register float ss;
1350  d=(1L<<max_height)-1-n;
1351  if (d>8) panic(bad_specs+1); /* too many vertices */
1352  if (d<0) nverts=0;
1353  else {
1354    nn[0]=nn[1]=1;
1355    for (k=2;k<=d;k++) nn[k]=0;
1356    for (j=2;j<=max_height;j++) {
1357      for (k=d;k;k--) {
1358        for (ss=0.0,i=k;i>=0;i--) ss+=((float)nn[i])*((float)nn[k-i]);
1359        if (ss>MAX_NNN) panic(very_bad_specs+1); /* way too big */
1360        for (s=0,i=k;i>=0;i--) s+=nn[i]*nn[k-i]; /* overflow impossible */
1361        nn[k]=s;
1362      }
1363      i=(1L<<j)-1;
1364      if (i<=d) nn[i]++; /* add $z^{1-2^j}$ */
1365    }
1366    nverts=nn[d];
1367  }
1368}
1369
1370@ We generate the trees in lexicographic order of their Polish prefix
1371notation, encoded in binary notation as $x_0x_1\ldots x_{2n}$, using `1'
1372for an internal node and `0' for a leaf. For example, the five
1373trees when $n=3$ are
1374$$1010100,\quad 1011000,\quad 1100100,\quad 1101000,\quad 1110000,$$
1375in lexicographic order. The algorithm for lexicographic generation maintains
1376three auxiliary arrays $l_j$, $y_j$, and $\sigma_j$, where
1377$$\sigma_j\;=\;n-j+\sum_{i=0}^{j-1}x_i\;=\;-1+\sum_{i=j}^{2n}(1-x_i)$$
1378is one less than the number of 0's (leaves) in $(x_j,\ldots,x_{2n})$.
1379The values of $l_j$ and $y_j$ are harder
1380to describe formally; $l_j$ is $2^{h-l}$ when $h=|max_height|$ and when
1381$x_j$ represents a node at level~$l$ of the tree, based on the values
1382of $(x_0,\ldots,x_{j-1})$. The value of $y_j$ is a binary encoding of
1383tree levels in which an internal node has not yet received a right child;
1384$y_j$ is also the maximum number of future leaves that can be produced by
1385previously specified internal nodes without exceeding the maximum height.
1386The number of 1-bits in $y_j$ is the minimum number of future leaves,
1387based on previous specifications.
1388
1389Therefore if $\sigma_j>y_j$, $x_j$ is forced to be~1. If $l_j=1$,
1390$x_j$ is forced to be~0. If the number of 1-bits of $y_j$ is equal
1391to $\sigma_j$, $x_j$ is forced to be~0. Otherwise $x_j$ can be
1392either 0 or~1, and it will be possible to complete the partial
1393solution $x_0\ldots x_j$ to a full Polish prefix code $x_0\ldots x_{2n}$.
1394
1395For example, here are the arrays for one of the binary trees
1396that is generated when $n=h=3$:
1397$$\vcenter{\halign{$\hfil#$\quad=&&\quad#\cr
1398j       &0&1&2&3&4&5&6\cr
1399l_j     &8&4&2&2&1&1&4\cr
1400y_j     &0&4&6&4&5&4&0\cr
1401\sigma_j&3&3&3&2&2&1&0\cr
1402x_j     &1&1&0&1&0&0&0\cr}}$$
1403If $x_j=1$ and $j<2n$, we have $l_{j+1}=l_j/2$, $y_{j+1}=y_j+l_{j+1}$,
1404and $\sigma_{j+1}=\sigma_j$. If $x_j=0$ and $j<2n$, we have $l_{j+1}=
14052^t$, $y_{j+1}=y_j-2^t$, and $\sigma_{j+1}=\sigma_j-1$, where $2^t$ is the
1406least power of~2 in the binary representation of~$y_j$. It is not difficult to
1407prove by induction that $\sigma_j<y_j+l_j$, assuming that $n<2^h$.
1408
1409@<Name the trees and create the arcs or edges@>=
1410{@+register long *xtab,*ytab,*ltab,*stab;
1411  @<Initialize |xtab|, |ytab|, |ltab|, and |stab|; also set |d=2n|@>;
1412  v=new_graph->vertices;
1413  if (ltab[0]>n) {
1414    k=0;@+xtab[0]=n?1:0;
1415    while (1) {
1416      @<Complete the partial tree $x_0\ldots x_k$@>;
1417      @<Assign a Polish prefix code name to vertex~|v|@>;
1418      @<Create arcs or edges from |v| to previous trees@>;
1419      v++;
1420      @<Advance to the next partial tree $x_0\ldots x_k$, where |k| is
1421       as large as possible; |goto last| if there are no more solutions@>;
1422    }
1423  }
1424}
1425last:@+if (v!=new_graph->vertices+new_graph->n)
1426  panic(impossible); /* can't happen */
1427gb_free(working_storage);
1428
1429@ @<Initialize |xtab|, |ytab|, |ltab|, and |stab|...@>=
1430xtab=gb_typed_alloc(8*n+4,long,working_storage);
1431if (gb_trouble_code) { /* no room for |xtab| */
1432  gb_recycle(new_graph);@+panic(no_room+2);@+}
1433d=n+n;
1434ytab=xtab+(d+1);
1435ltab=ytab+(d+1);
1436stab=ltab+(d+1);
1437ltab[0]=1L<<max_height;
1438stab[0]=n; /* |ytab[0]=0| */
1439
1440@ @<Complete the partial tree...@>=
1441for (j=k+1;j<=d;j++) {
1442  if (xtab[j-1]) {
1443    ltab[j]=ltab[j-1]>>1;
1444    ytab[j]=ytab[j-1]+ltab[j];
1445    stab[j]=stab[j-1];
1446  }@+else {
1447    ytab[j]=ytab[j-1]&(ytab[j-1]-1); /* remove least significant 1-bit */
1448    ltab[j]=ytab[j-1]-ytab[j];
1449    stab[j]=stab[j-1]-1;
1450  }
1451  if (stab[j]<=ytab[j]) xtab[j]=0;
1452  else xtab[j]=1;  /* this is the lexicographically smallest completion */
1453}
1454
1455@ As in previous routines, we seek the largest $k$ such that $x_k$ can
1456be increased without violating the necessary and sufficient conditions
1457stated earlier.
1458
1459@<Advance to the next partial tree...@>=
1460for (k=d-1;;k--) {
1461  if (k<=0) goto last; /* this happens only when |n<=1| */
1462  if (xtab[k]) break; /* find rightmost 1 */
1463}
1464for (k--;;k--) {
1465  if (xtab[k]==0 && ltab[k]>1) break;
1466  if (k==0) goto last;
1467}
1468xtab[k]++;
1469
1470@ In the |name| field, we encode internal nodes of the binary tree by
1471`\..' and leaves by `\.x'. Thus the five trees shown above in binary
1472code will be named
1473$$\.{.x.x.xx},\quad \.{.x..xxx},\quad \.{..xx.xx},\quad \.{..x.xxx},\quad
1474\.{...xxxx},$$
1475respectively.
1476
1477@<Assign a Polish prefix...@>=
1478{@+register char *p=buffer; /* string pointer */
1479  for (k=0;k<=d;k++,p++) *p=(xtab[k]? '.': 'x');
1480  v->name=gb_save_string(buffer);
1481  hash_in(v); /* enter |v->name| into the hash table
1482                  (via utility fields |u,v|) */
1483}
1484
1485@ Since we are generating the trees in lexicographic order of their
1486Polish prefix notation, it is relatively easy to find all pairs of trees that
1487are adjacent via one application of the associative law: We simply
1488replace a substring of the form $\..\..\alpha\beta$ by
1489$\..\alpha\..\beta$, when $\alpha$ and $\beta$ are Polish prefix
1490strings. The result comes earlier in lexicographic order, so it will
1491be an existing vertex unless it violates the |max_height| restriction.
1492
1493@<Create arcs or edges from |v| to previous trees@>=
1494for (j=0;j<d;j++)
1495  if (xtab[j]==1 && xtab[j+1]==1) {
1496    for (i=j+1,s=0;s>=0;s+=(xtab[i+1]<<1)-1,i++) xtab[i]=xtab[i+1];
1497    xtab[i]=1;
1498    {@+register char *p=buffer; /* string pointer */
1499      register Vertex *u;
1500      for (k=0;k<=d;k++,p++) *p=(xtab[k]? '.': 'x');
1501      u=hash_out(buffer);
1502      if (u) {
1503        if (directed) gb_new_arc(v,u,1L);
1504        else gb_new_edge(v,u,1L);
1505      }
1506    }
1507    for (i--;i>j;i--) xtab[i+1]=xtab[i]; /* restore |xtab| */
1508    xtab[i+1]=1;
1509  }
1510
1511@* Complementing and copying. We have seen how to create a wide
1512variety of basic graphs with the |board|, |simplex|, |subsets|,
1513|perms|, |parts|, and |binary| procedures. The remaining routines
1514of {\sc GB\_\,BASIC} are somewhat different. They transform existing
1515graphs into new ones, thereby presenting us with an almost
1516mind-boggling array of further possibilities.
1517
1518The first of these transformations is perhaps the simplest. It
1519complements a given graph, making vertices adjacent if and only if
1520they were previously nonadjacent. More precisely, the subroutine call
1521|complement(g,copy,self,directed)| returns a graph with the
1522same vertices as |g|, but with complemented arcs.
1523If |self!=0|, the new graph will have a self-loop from a vertex |v| to itself
1524when the original graph did not; if |self=0|, the new graph will
1525have no self-loops. If |directed!=0|, the new graph will have
1526an arc from |u| to |v| when the original graph did not; if |directed=0|,
1527the new graph will be undirected, and it will have an edge between |u|
1528and~|v| when the original graph did not. In the latter case, the original
1529graph should also be undirected (that is, its arcs should come in pairs,
1530as described in the |gb_new_edge| routine of {\sc GB\_\,GRAPH}).
1531
1532If |copy!=0|, a double complement will actually be done. This means that
1533the new graph will essentially be a copy of the old, except that
1534duplicate arcs (and possibly self-loops) will be removed. Regardless of
1535the value of |copy|, information that might have been present in the utility
1536fields will not be copied, and arc lengths will all be set to~1.
1537
1538One possibly useful feature of the graphs returned by |complement| is
1539worth noting. The vertices adjacent to~|v|, namely the list
1540$$\hbox{|v->arcs->tip|,\quad |v->arcs->next->tip|,\quad
1541 |v->arcs->next->next->tip|,\quad \dots\thinspace,}$$
1542will be in strictly decreasing order (except in the case of an
1543undirected self-loop, when |v| itself will appear twice in succession).
1544
1545@ @<Basic subroutines@>=
1546Graph *complement(g,copy,self,directed)
1547  Graph *g; /* graph to be complemented */
1548  long copy; /* should we double-complement? */
1549  long self; /* should we produce self-loops? */
1550  long directed; /* should the graph be directed? */
1551{@+@<Vanilla local variables@>@;@+@q{@>@t}\6{@>@q}@>
1552  register long n;
1553  register Vertex *u;
1554  register siz_t delta; /* difference in memory addresses */
1555  if (g==NULL) panic(missing_operand); /* where's |g|? */
1556  @<Set up a graph with the vertices of |g|@>;
1557  sprintf(buffer,",%d,%d,%d)",copy?1:0,self?1:0,directed?1:0);
1558  make_compound_id(new_graph,"complement(",g,buffer);
1559  @<Insert complementary arcs or edges@>;
1560  if (gb_trouble_code) {
1561    gb_recycle(new_graph);
1562    panic(alloc_fault);
1563     /* worse luck, we ran out of memory somewhere back there */
1564  }
1565  return new_graph;
1566}
1567
1568@ In several of the following routines, it is efficient to circumvent
1569\CEE/'s normal rules for pointer arithmetic, and to use the
1570fact that the vertices of a graph being copied are a constant distance away
1571in memory from the vertices of its clone.
1572
1573@d vert_offset(v,delta) ((Vertex*)(((siz_t)v)+delta))
1574@^pointer hacks@>
1575
1576@<Set up a graph with the vertices of |g|@>=
1577n=g->n;
1578new_graph=gb_new_graph(n);
1579if (new_graph==NULL)
1580  panic(no_room); /* out of memory before we're even started */
1581delta=((siz_t)(new_graph->vertices))-((siz_t)(g->vertices));
1582for (u=new_graph->vertices,v=g->vertices;v<g->vertices+n;u++,v++)
1583  u->name=gb_save_string(v->name);
1584
1585@ A temporary utility field in the new graph is used to remember which
1586vertices are adjacent to a given vertex in the old one. We stamp the |tmp|
1587field of~|v| with a pointer to~|u| when there's an arc from |u| to~|v|.
1588
1589@d tmp u.V /* utility field |u| for temporary use as a vertex pointer */
1590
1591@<Insert comp...@>=
1592for (v=g->vertices;v<g->vertices+n;v++) {@+register Vertex *vv;
1593  u=vert_offset(v,delta);
1594           /* vertex in |new_graph| corresponding to |v| in |g| */
1595  {@+register Arc *a;
1596    for (a=v->arcs;a;a=a->next) vert_offset(a->tip,delta)->tmp=u;
1597  }
1598  if (directed) {
1599    for (vv=new_graph->vertices;vv<new_graph->vertices+n;vv++)
1600      if ((vv->tmp==u && copy) || (vv->tmp!=u && !copy))
1601        if (vv!=u || self) gb_new_arc(u,vv,1L);
1602  }@+else {
1603    for (vv=(self?u:u+1);vv<new_graph->vertices+n;vv++)
1604      if ((vv->tmp==u && copy) || (vv->tmp!=u && !copy))
1605        gb_new_edge(u,vv,1L);
1606  }
1607}
1608for (v=new_graph->vertices;v<new_graph->vertices+n;v++) v->tmp=NULL;
1609
1610@* Graph union and intersection. Another simple way to get new graphs
1611from old ones is to take the union or intersection of their sets of arcs. The
1612subroutine call |gunion(g,gg,multi,directed)| produces a graph
1613with the vertices and arcs of |g| together with the
1614arcs of another graph~|gg|. The subroutine call |intersection(g,gg,multi,
1615directed)| produces a graph with the vertices of |g| but with only the
1616arcs that appear in both |g| and |gg|. In both cases we assume
1617that |gg| has the same vertices as |g|, in the sense that vertices
1618in the same relative position from the beginning of the vertex array
1619are considered identical. If the actual number of vertices in |gg| exceeds
1620the number in |g|, the extra vertices and all arcs touching them in~|gg| are
1621suppressed.
1622
1623The input graphs are assumed to be undirected, unless the |directed|
1624parameter is nonzero. Peculiar results might occur if you mix directed
1625and undirected graphs, but the subroutines will not ``crash''
1626when they are asked to produce undirected output from directed input.
1627
1628If |multi| is nonzero, the new graph may have multiple edges: Suppose
1629there are $k_1$ arcs from $u$ to $v$ in |g| and $k_2$ such arcs in |gg|. Then
1630there will be $k_1+k_2$ in the union and $\min(k_1,k_2)$ in the
1631intersection when |multi!=0|, but at most
1632one in the union or intersection when |multi=0|.
1633
1634The lengths of arcs are copied to the union graph when |multi!=0|;
1635the minimum length of multiple arcs is retained in the union when |multi=0|.
1636
1637The lengths of arcs in the intersection graph are a bit trickier.
1638If multiple arcs occur in |g|, their minimum length, |l|, is computed. Then
1639we compute the maximum of |l| and the lengths of corresponding arcs
1640in |gg|. If |multi=0|, only the minimum of those maxima will survive.
1641
1642@ @<Basic subroutines@>=
1643Graph *gunion(g,gg,multi,directed)
1644  Graph *g,*gg; /* graphs to be united */
1645  long multi; /* should we reproduce multiple arcs? */
1646  long directed; /* should the graph be directed? */
1647{@+@<Vanilla local variables@>@;@+@q{@>@t}\6{@>@q}@>
1648  register long n;
1649  register Vertex *u;
1650  register siz_t delta,ddelta; /* differences in memory addresses */
1651  if (g==NULL || gg==NULL) panic(missing_operand);
1652    /* where are |g| and |gg|? */
1653  @<Set up a graph with the vertices of |g|@>;
1654  sprintf(buffer,",%d,%d)",multi?1:0,directed?1:0);
1655  make_double_compound_id(new_graph,"gunion(",g,",",gg,buffer);
1656  ddelta=((siz_t)(new_graph->vertices))-
1657             ((siz_t)(gg->vertices));
1658  @<Insert arcs or edges present in either |g| or |gg|@>;
1659  if (gb_trouble_code) {
1660    gb_recycle(new_graph);
1661    panic(alloc_fault); /* uh oh, we ran out of memory somewhere back there */
1662  }
1663  return new_graph;
1664}
1665
1666@ @<Insert arcs or edges present in either |g| or |gg|@>=
1667for (v=g->vertices;v<g->vertices+n;v++) {
1668  register Arc *a;
1669  register Vertex *vv=vert_offset(v,delta);
1670      /* vertex in |new_graph| corresponding to |v| in |g| */
1671  register Vertex *vvv=vert_offset(vv,-ddelta);
1672      /* vertex in |gg| corresponding to |v| in  |g| */
1673  for (a=v->arcs;a;a=a->next) {
1674    u=vert_offset(a->tip,delta);
1675    @<Insert a union arc or edge from |vv| to |u|, if appropriate@>;
1676  }
1677  if (vvv<gg->vertices+gg->n) for (a=vvv->arcs;a;a=a->next) {
1678    u=vert_offset(a->tip,ddelta);
1679    if (u<new_graph->vertices+n)
1680      @<Insert a union arc or edge from |vv| to |u|, if appropriate@>;
1681  }
1682}
1683for (v=new_graph->vertices;v<new_graph->vertices+n;v++)
1684  v->tmp=NULL,v->tlen=NULL;
1685
1686@ We use the |tmp| trick of |complement| to remember which arcs have
1687already been recorded from |u|, and we extend it so that we can maintain
1688minimum lengths. Namely, |uu->tmp| will equal |u| if and only
1689if we have already seen an arc from |u| to |uu|; and if so, |uu->tlen|
1690will be one such arc. In the undirected case, |uu->tlen| will point to the
1691first arc of an edge pair that touches~|u|.
1692
1693The only thing slightly nontrivial here is the way we keep undirected
1694edges grouped in pairs. We generate a new edge from |vv| to |u| only
1695if |vv<=u|, and if equality holds we advance~|a| so that we don't
1696see the self-loop in both directions. Similar logic will be repeated
1697in many of the programs below.
1698
1699@d tlen z.A /* utility field |z| regarded as a pointer to an arc */
1700
1701@<Insert a union arc or edge from |vv| to |u|, if appropriate@>=
1702{@+register Arc *b;
1703  if (directed) {
1704    if (multi || u->tmp!=vv) gb_new_arc(vv,u,a->len);
1705    else {
1706      b=u->tlen;
1707      if (a->len<b->len) b->len=a->len;
1708    }
1709    u->tmp=vv; /* remember that we've seen this */
1710    u->tlen=vv->arcs;
1711  }@+else if (u>=vv) {
1712    if (multi || u->tmp!=vv) gb_new_edge(vv,u,a->len);
1713    else {
1714      b=u->tlen;
1715      if (a->len<b->len) b->len=(b+1)->len=a->len;
1716    }
1717    u->tmp=vv;
1718    u->tlen=vv->arcs;
1719    if (u==vv && a->next==a+1) a++; /* bypass second half of self-loop */
1720  }
1721}
1722
1723@ @<Basic subroutines@>=
1724Graph *intersection(g,gg,multi,directed)
1725  Graph *g,*gg; /* graphs to be intersected */
1726  long multi; /* should we reproduce multiple arcs? */
1727  long directed; /* should the graph be directed? */
1728{@+@<Vanilla local variables@>@;@+@q{@>@t}\6{@>@q}@>
1729  register long n;
1730  register Vertex *u;
1731  register siz_t delta,ddelta; /* differences in memory addresses */
1732  if (g==NULL || gg==NULL) panic(missing_operand); /* where are |g| and |gg|? */
1733  @<Set up a graph with the vertices of |g|@>;
1734  sprintf(buffer,",%d,%d)",multi?1:0,directed?1:0);
1735  make_double_compound_id(new_graph,"intersection(",g,",",gg,buffer);
1736  ddelta=((siz_t)(new_graph->vertices))-
1737               ((siz_t)(gg->vertices));
1738  @<Insert arcs or edges present in both |g| and |gg|@>;
1739  if (gb_trouble_code) {
1740    gb_recycle(new_graph);
1741    panic(alloc_fault); /* whoops, we ran out of memory somewhere back there */
1742  }
1743  return new_graph;
1744}
1745
1746@ Two more temporary utility fields are needed here.
1747
1748@d mult v.I /* utility field |v|, counts multiplicity of arcs */
1749@d minlen w.I /* utility field |w|, records the smallest length */
1750
1751@<Insert arcs or edges present in both |g| and |gg|@>=
1752for (v=g->vertices;v<g->vertices+n;v++) {@+register Arc *a;
1753  register Vertex *vv=vert_offset(v,delta);
1754      /* vertex in |new_graph| corresponding to |v| in |g| */
1755  register Vertex *vvv=vert_offset(vv,-ddelta);
1756      /* vertex in |gg| corresponding to |v| in  |g| */
1757  if (vvv>=gg->vertices+gg->n) continue;
1758  @<Take note of all arcs from |v|@>;
1759  for (a=vvv->arcs;a;a=a->next) {
1760    u=vert_offset(a->tip,ddelta);
1761    if (u>=new_graph->vertices+n) continue;
1762    if (u->tmp==vv) {@+long l=u->minlen;
1763      if (a->len>l) l=a->len; /* maximum */
1764      if (u->mult<0) @<Update minimum of multiple maxima@>@;
1765      else @<Generate a new arc or edge for the intersection,
1766                and reduce the multiplicity@>;
1767    }
1768  }
1769}
1770@<Clear out the temporary utility fields@>;
1771
1772@ @<Generate a new arc or edge for the intersection...@>=
1773{
1774  if (directed) gb_new_arc(vv,u,l);
1775  else {
1776    if (vv<=u) gb_new_edge(vv,u,l);
1777    if (vv==u && a->next==a+1) a++; /* skip second half of self-loop */
1778  }
1779  if (!multi) {
1780    u->tlen=vv->arcs;
1781    u->mult=-1;
1782  }@+else if (u->mult==0) u->tmp=NULL;
1783  else u->mult--;
1784}
1785
1786@ We get here if and only if |multi=0| and |gg|~has more than one arc from |vv|
1787to~|u| and |g|~has at least one arc from |vv| to~|u|.
1788
1789@<Update minimum of multiple maxima@>=
1790{@+register Arc *b=u->tlen; /* previous arc or edge from |vv| to |u| */
1791  if (l<b->len) {
1792    b->len=l;
1793    if (!directed) (b+1)->len=l;
1794  }
1795}
1796
1797@ @<Take note of all arcs from |v|@>=
1798for (a=v->arcs;a;a=a->next) {
1799  u=vert_offset(a->tip,delta);
1800  if (u->tmp==vv) {
1801    u->mult++;
1802    if (a->len<u->minlen) u->minlen=a->len;
1803  }@+else u->tmp=vv, u->mult=0, u->minlen=a->len;
1804  if (u==vv && !directed && a->next==a+1) a++;
1805           /* skip second half of self-loop */
1806}
1807
1808@ @<Clear out the temporary utility fields@>=
1809for (v=new_graph->vertices;v<new_graph->vertices+n;v++) {
1810  v->tmp=NULL;
1811  v->tlen=NULL;
1812  v->mult=0;
1813  v->minlen=0;
1814}
1815
1816@* Line graphs. The next operation in {\sc GB\_\,BASIC}'s repertoire constructs
1817the so-called line graph of a given graph~$g$. The subroutine that does
1818this is invoked by calling `|lines(g,directed)|'.
1819
1820If |directed=0|, the line graph has one vertex for each edge of~|g|;
1821two vertices are adjacent if and only if the corresponding edges
1822have a common vertex.
1823
1824If |directed!=0|, the line graph has one vertex for each arc of~|g|;
1825there is an arc from vertex |u| to vertex |v| if and only if the
1826arc corresponding to~|u| ends at the vertex that begins the arc
1827corresponding to~|v|.
1828
1829All arcs of the line graph will have length~1.
1830
1831Utility fields |u.V| and |v.V| of each vertex in the line graph will point to
1832the vertices of |g| that define the corresponding arc or edge, and |w.A| will
1833point to the arc from |u.V| to |v.V| in~|g|. In the undirected case we will
1834have |u.V<=v.V|.
1835
1836@<Basic subroutines@>=
1837Graph *lines(g,directed)
1838  Graph *g; /* graph whose lines will become vertices */
1839  long directed; /* should the graph be directed? */
1840{@+@<Vanilla local variables@>@;@+@q{@>@t}\6{@>@q}@>
1841  register long m; /* the number of lines */
1842  register Vertex *u;
1843  if (g==NULL) panic(missing_operand); /* where is |g|? */
1844  @<Set up a graph whose vertices are the lines of |g|@>;
1845  if (directed) @<Insert arcs of a directed line graph@>@;
1846  else @<Insert edges of an undirected line graph@>;
1847  @<Restore |g| to its pristine original condition@>;
1848  if (gb_trouble_code) {
1849    gb_recycle(new_graph);
1850    panic(alloc_fault); /* (sigh) we ran out of memory somewhere back there */
1851  }
1852  return new_graph;
1853near_panic:@<Recover from potential disaster due to bad data@>;
1854}
1855
1856@ We want to add a data structure to |g| so that the line graph can be
1857built efficiently. But we also want to preserve |g| so that it
1858exhibits no traces of occupation when |lines| has finished its
1859work.  To do this, we will move utility field~|v->z| temporarily into
1860a utility field~|u->z| of the line graph, where |u| is the first
1861vertex having |u->u.V==v|, whenever such a |u| exists. Then we'll
1862set |v->map=u|. We will then be able to find |u| when |v|
1863is given, and we'll be able to cover our tracks later.
1864
1865In the undirected case, further structure is needed. We will temporarily
1866change the |tip| field in the second arc of each edge pair so that
1867it points to the line-graph vertex that points to the first arc of the pair.
1868
1869The |util_types| field of the graph does not indicate the fact that utility
1870fields |u.V|, |v.V|, and |w.A| of each vertex will be set, because those
1871utility fields are pointers from the new graph to the original graph.
1872The |save_graph| procedure does not deal with pointers between graphs.
1873
1874@d map z.V /* the |z| field treated as a vertex pointer */
1875
1876@<Restore |g| to its pristine original condition@>=
1877for (u=new_graph->vertices,v=NULL;u<new_graph->vertices+m;u++) {
1878  if (u->u.V!=v) {
1879    v=u->u.V; /* original vertex of |g| */
1880    v->map=u->map; /* restore original value of |v->z| */
1881    u->map=NULL;
1882  }
1883  if (!directed) ((u->w.A)+1)->tip=v;
1884}
1885
1886@ Special care must be taken to avoid chaos when the user is trying to
1887construct the undirected line graph of a directed graph.  Otherwise we
1888might trash the memory, or we might leave the original graph in a garbled state
1889with pointers leading into supposedly free space.
1890
1891@<Set up a graph whose vertices are the lines of |g|@>=
1892m=(directed? g->m: (g->m)/2);
1893new_graph=gb_new_graph(m);
1894if (new_graph==NULL)
1895  panic(no_room); /* out of memory before we're even started */
1896make_compound_id(new_graph,"lines(",g,directed? ",1)": ",0)");
1897u=new_graph->vertices;
1898for (v=g->vertices+g->n-1;v>=g->vertices;v--) {@+register Arc *a;
1899  register long mapped=0; /* has |v->map| been set? */
1900  for (a=v->arcs;a;a=a->next) {@+register Vertex *vv=a->tip;
1901    if (!directed) {
1902      if (vv<v) continue;
1903      if (vv>=g->vertices+g->n) goto near_panic;
1904            /* original graph not undirected */
1905    }
1906    @<Make |u| a vertex representing the arc |a| from |v| to |vv|@>;
1907    if (!mapped) {
1908      u->map=v->map;  /* |z.V=map| incorporates all bits of utility field |z|,
1909                           whatever its type */
1910      v->map=u;
1911      mapped=1;
1912    }
1913    u++;
1914  }
1915}
1916if (u!=new_graph->vertices+m) goto near_panic;
1917
1918@ @<Recover...@>=
1919m=u-new_graph->vertices;
1920@<Restore |g| to its pristine...@>;
1921gb_recycle(new_graph);
1922panic(invalid_operand);
1923 /* |g| did not obey the conventions for an undirected graph */
1924
1925@ The vertex names in the line graph are pairs of original vertex names,
1926separated by `\.{--}' when undirected, `\.{->}' when directed. If either
1927of the original names is horrendously long, the villainous Procrustes
1928chops it off arbitrarily so that it fills at most half of the name buffer.
1929
1930@<Make |u| a vertex representing the arc |a| from |v| to |vv|@>=
1931u->u.V=v;
1932u->v.V=vv;
1933u->w.A=a;
1934if (!directed) {
1935  if (u>=new_graph->vertices+m || (a+1)->tip!=v) goto near_panic;
1936  if (v==vv && a->next==a+1) a++; /* skip second half of self-loop */
1937  else (a+1)->tip=u;
1938}
1939sprintf(buffer,"%.*s-%c%.*s",(BUF_SIZE-3)/2,v->name,@|
1940   directed? '>': '-',BUF_SIZE/2-1,vv->name);
1941u->name=gb_save_string(buffer);
1942
1943@ @<Insert arcs of a directed line graph@>=
1944for (u=new_graph->vertices;u<new_graph->vertices+m;u++) {
1945  v=u->v.V;
1946  if (v->arcs) { /* |v->map| has been set up */
1947    v=v->map;
1948    do@+{gb_new_arc(u,v,1L);
1949    v++;
1950    }@+while (v->u.V==u->v.V);
1951  }
1952}
1953
1954@ An undirected line graph will contain no self-loops. It contains
1955multiple edges only if the original graph did; in that case, there
1956are two edges joining a line to each of its parallel mates, because
1957each mate hits both of its endpoints.
1958
1959The details of this section deserve careful study.  We use the
1960fact that the first vertices of the lines occur in nonincreasing order.
1961
1962@<Insert edges of an undirected line graph@>=
1963for (u=new_graph->vertices;u<new_graph->vertices+m;u++) {@+register Vertex *vv;
1964  register Arc *a;@+register long mapped=0;
1965  v=u->u.V; /* we look first for prior lines that touch the first vertex */
1966  for (vv=v->map;vv<u;vv++) gb_new_edge(u,vv,1L);
1967  v=u->v.V; /* then we look for prior lines that touch the other one */
1968  for (a=v->arcs;a;a=a->next) {
1969    vv=a->tip;
1970    if (vv<u && vv>=new_graph->vertices) gb_new_edge(u,vv,1L);
1971    else if (vv>=v && vv<g->vertices+g->n) mapped=1;
1972  }
1973  if (mapped && v>u->u.V)
1974    for (vv=v->map;vv->u.V==v;vv++) gb_new_edge(u,vv,1L);
1975}
1976
1977@* Graph products. Three ways have traditionally been used to define the
1978product of two graphs. In all three cases the vertices of the product graph
1979are ordered pairs $(v,v')$, where $v$ and $v'$ are vertices of the original
1980graphs; the difference occurs in the definition of arcs. Suppose $g$ has
1981$m$ arcs and $n$ vertices, while $g'$ has $m'$ arcs and $n'$ vertices. The
1982{\sl cartesian product\/} of $g$ and~$g'$ has $mn'+m'n$ arcs, namely from
1983$(u,u')$ to $(v,u')$ whenever there's an arc from $u$ to $v$ in~$g$, and from
1984$(u,u')$ to $(u,v')$ whenever there's an arc from $u'$ to $v'$ in~$g'$.
1985The {\sl direct product\/} has $mm'$ arcs, namely from $(u,u')$ to
1986$(v,v')$ in the same circumstances. The {\sl strong product\/}
1987has both the arcs of the cartesian product and the direct product.
1988
1989Notice that an undirected graph with $m$ edges has $2m$ arcs. Thus the
1990number of edges in the direct product of two undirected graphs is
1991twice the product of the number of edges in the individual graphs.
1992A self-loop in~$g$ will combine with an edge in~$g'$ to make
1993two parallel edges in the direct product.
1994
1995The subroutine call |product(g,gg,type,directed)| produces the product
1996graph of one of these three types, where |type=0| for cartesian product,
1997|type=1| for direct product, and |type=2| for strong product.
1998The length of an arc in the cartesian product is copied from the length
1999of the original arc that it replicates; the length of an arc in the direct
2000product is the minimum of the two arc lengths that induce it. If |directed=0|,
2001the product graph will be an undirected graph with edges consisting
2002of consecutive arc pairs according to the standard GraphBase conventions,
2003and the input graphs should adhere to  the same conventions.
2004
2005@(gb_basic.h@>=
2006#define cartesian 0
2007#define direct 1
2008#define strong 2
2009
2010@ @<Basic subroutines@>=
2011Graph *product(g,gg,type,directed)
2012  Graph *g,*gg; /* graphs to be multiplied */
2013  long type; /* |cartesian|, |direct|, or |strong| */
2014  long directed; /* should the graph be directed? */
2015{@+@<Vanilla local variables@>@;@+@q{@>@t}\6{@>@q}@>
2016  register Vertex *u,*vv;
2017  register long n; /* the number of vertices in the product graph */
2018  if (g==NULL || gg==NULL) panic(missing_operand); /* where are |g| and |gg|? */
2019  @<Set up a graph with ordered pairs of vertices@>;
2020  if ((type&1)==0) @<Insert arcs or edges for cartesian product@>;
2021  if (type) @<Insert arcs or edges for direct product@>;
2022  if (gb_trouble_code) {
2023    gb_recycle(new_graph);
2024    panic(alloc_fault);
2025      /* @@?`$*$\#!\&, we ran out of memory somewhere back there */
2026  }
2027  return new_graph;
2028}
2029
2030@ We must be constantly on guard against running out of memory, especially
2031when multiplying information.
2032
2033The vertex names in the product are pairs of original vertex names separated
2034by commas. Thus, for example, if you cross an |econ| graph with a |roget|
2035graph, you can get vertices like |"Financial services, mediocrity"|.
2036
2037@<Set up a graph with ordered pairs of vertices@>=
2038{@+float test_product=((float)(g->n))*((float)(gg->n));
2039  if (test_product>MAX_NNN) panic(very_bad_specs); /* way too many vertices */
2040}
2041n=(g->n)*(gg->n);
2042new_graph=gb_new_graph(n);
2043if (new_graph==NULL)
2044  panic(no_room); /* out of memory before we're even started */
2045for (u=new_graph->vertices,v=g->vertices,vv=gg->vertices;@|
2046    u<new_graph->vertices+n;u++) {
2047  sprintf(buffer,"%.*s,%.*s",BUF_SIZE/2-1,v->name,(BUF_SIZE-1)/2,vv->name);
2048  u->name=gb_save_string(buffer);
2049  if (++vv==gg->vertices+gg->n) vv=gg->vertices,v++; /* ``carry'' */
2050}
2051sprintf(buffer,",%d,%d)",(type?2:0)-(int)(type&1),directed?1:0);
2052make_double_compound_id(new_graph,"product(",g,",",gg,buffer);
2053
2054@ @<Insert arcs or edges for cartesian product@>=
2055{@+register Vertex *uu,*uuu;
2056  register Arc *a;
2057  register siz_t delta; /* difference in memory addresses */
2058  delta=((siz_t)(new_graph->vertices))-((siz_t)(gg->vertices));
2059  for (u=gg->vertices;u<gg->vertices+gg->n;u++)
2060    for (a=u->arcs;a;a=a->next) {
2061      v=a->tip;
2062      if (!directed) {
2063        if (u>v) continue;
2064        if (u==v && a->next==a+1) a++; /* skip second half of self-loop */
2065      }
2066      for (uu=vert_offset(u,delta),vv=vert_offset(v,delta);@|
2067             uu<new_graph->vertices+n;uu+=gg->n,vv+=gg->n)
2068        if (directed) gb_new_arc(uu,vv,a->len);
2069        else gb_new_edge(uu,vv,a->len);
2070    }
2071  @<Insert arcs or edges for first component of cartesian product@>;
2072}
2073
2074@ @<Insert arcs or edges for first component...@>=
2075for (u=g->vertices,uu=new_graph->vertices;uu<new_graph->vertices+n;
2076          u++,uu+=gg->n)
2077  for (a=u->arcs;a;a=a->next) {
2078    v=a->tip;
2079    if (!directed) {
2080      if (u>v) continue;
2081      if (u==v && a->next==a+1) a++; /* skip second half of self-loop */
2082    }
2083    vv=new_graph->vertices+((gg->n)*(v-g->vertices));
2084    for (uuu=uu;uuu<uu+gg->n;uuu++,vv++)
2085      if (directed) gb_new_arc(uuu,vv,a->len);
2086      else gb_new_edge(uuu,vv,a->len);
2087  }
2088
2089@ @<Insert arcs or edges for direct product@>=
2090{@+Vertex *uu;@+Arc *a;
2091  siz_t delta0=
2092   ((siz_t)(new_graph->vertices))-((siz_t)(gg->vertices));
2093  siz_t del=(gg->n)*sizeof(Vertex);
2094  register siz_t delta,ddelta;
2095  for (uu=g->vertices,delta=delta0;uu<g->vertices+g->n;uu++,delta+=del)
2096    for (a=uu->arcs;a;a=a->next) {
2097      vv=a->tip;
2098      if (!directed) {
2099        if (uu>vv) continue;
2100        if (uu==vv && a->next==a+1) a++; /* skip second half of self-loop */
2101      }
2102      ddelta=delta0+del*(vv-g->vertices);
2103      for (u=gg->vertices;u<gg->vertices+gg->n;u++) {@+register Arc *aa;
2104        for (aa=u->arcs;aa;aa=aa->next) {@+long length=a->len;
2105          if (length>aa->len) length=aa->len;
2106          v=aa->tip;
2107          if (directed)
2108            gb_new_arc(vert_offset(u,delta),vert_offset(v,ddelta),length);
2109          else gb_new_edge(vert_offset(u,delta),
2110                           vert_offset(v,ddelta),length);
2111        }
2112      }
2113    }
2114}
2115
2116@* Induced graphs. Another important way to transform a graph is to
2117remove, identify, or split some of its vertices. All of these
2118operations are performed by the |induced| routine, which users can
2119invoke by calling `|induced(g,description,self,multi,directed)|'.
2120
2121Each vertex |v| of |g| should first be assigned an ``induction code'' in
2122its field |v->ind|, which is actually utility field~|z|. The
2123induction code is 0~if |v| is to be eliminated; it is 1~if |v| is to be
2124retained; it is |k>1| if |v| is to be split into $k$ nonadjacent vertices
2125having the same neighbors as~|v| did; and it is |k<0| if |v| is to be
2126identified with all other vertices having the same value of~|k|.
2127
2128For example, suppose |g| is a circuit with vertices $\{0,1,\ldots,9\}$,
2129where |j| is adjacent to~|k| if and only if $k=(j\pm1)\bmod10$.
2130If we set
2131$$\vcenter{\halign{\hbox{\hfil#\hfil}\cr
2132|0->ind=0|,\quad |1->ind=5->ind=9->ind=-1|,\quad |2->ind=3->ind=-2|,\cr
2133|4->ind=6->ind=8->ind=1|,\quad and |7->ind=3|,\cr}}$$
2134the induced graph will have vertices $\{-1,-2,4,6,7,7',7'',8\}$.
2135The vertices adjacent to 6, say, will be $-1$ (formerly~5), 7, $7'$,
2136and~$7''$. The vertices adjacent to $-1$ will be those formerly
2137adjacent to 1, 5, or~9, namely $-2$ (formerly~2), 4, 6, and~8. The
2138vertices adjacent to $-2$ will be those formerly adjacent to 2 or~3,
2139namely $-1$ (formerly~1), $-2$ (formerly~3), $-2$ (formerly~2), and~4.
2140Duplicate edges will be discarded if |multi==0|, and self-loops will
2141be discarded if |self==0|.
2142
2143The total number of vertices in the induced graph will be the sum
2144of the positive |ind| fields plus the absolute value of the most
2145negative |ind| field. This rule implies, for example, that if at least
2146one vertex has |ind=-5|, the induced graph will always have a vertex $-4$,
2147even though no |ind| field has been set to~$-4$.
2148
2149The |description| parameter is a string that will appear as part of
2150the name of the induced graph; if |description=0|, this string will
2151be empty. In the latter case, users are encouraged to assign a suitable
2152name to the |id| field of the induced graph themselves, characterizing
2153the method by which the |ind| codes were set.
2154
2155If the |directed| parameter is zero, the input graph will be assumed to
2156be undirected, and the output graph will be undirected.
2157
2158When |multi=0|, the length of an arc that represents multiple arcs
2159will be the minimum of the multiple arc lengths.
2160
2161@d ind z.I
2162
2163@(gb_basic.h@>=
2164#define ind @[z.I /* utility field |z| when used to induce a graph */@]
2165
2166@ Here's a simple example: To get a complete bipartite graph with
2167parts of sizes |n1| and |n2|, we can start with a trivial two-point
2168graph and split its vertices into |n1| and |n2| parts.
2169
2170@<Applications...@>=
2171Graph *bi_complete(n1,n2,directed)
2172  unsigned long n1; /* size of first part */
2173  unsigned long n2; /* size of second part */
2174  long directed; /* should all arcs go from first part to second? */
2175{@+Graph *new_graph=board(2L,0L,0L,0L,1L,0L,directed);
2176  if (new_graph) {
2177    new_graph->vertices->ind=n1;
2178    (new_graph->vertices+1)->ind=n2;
2179    new_graph=induced(new_graph,NULL,0L,0L,directed);
2180    if (new_graph) {
2181      sprintf(new_graph->id,"bi_complete(%lu,%lu,%d)",@|
2182                              n1,n2,directed?1:0);
2183      mark_bipartite(new_graph,n1);
2184    }
2185  }
2186  return new_graph;
2187}
2188
2189@ The |induced| routine also provides a special feature not mentioned
2190above: If the |ind| field of any vertex |v| is |IND_GRAPH| or greater
2191(where |IND_GRAPH| is a large constant, much larger than the number
2192of vertices that would fit in computer memory), then utility field |v->subst|
2193should point to a graph. A copy of the vertices of
2194that graph will then be substituted for |v| in the induced graph.
2195
2196This feature extends the ordinary case when |v->ind>0|, which essentially
2197substitutes an empty graph for~|v|.
2198
2199If substitution is being used to replace all of $g$'s vertices
2200by disjoint copies of some other graph~$g'$,
2201the induced graph will be somewhat similar to
2202a product graph. But it will not be the same as any of the three
2203types of output produced by |product|, because the relation between
2204$g$ and $g'$ is not symmetrical. Assuming that no self-loops are
2205present, and that graphs $(g,g')$ have respectively $(m,m')$ arcs and
2206$(n,n')$ vertices, the result of substituting $g'$ for all
2207vertices of~$g$ has $m'n+mn'^2$ arcs.
2208
2209
2210@d IND_GRAPH 1000000000 /* when |ind| is a billion or more, */
2211@d subst y.G /* we'll look at the |subst| field */
2212
2213@(gb_basic.h@>=
2214#define IND_GRAPH 1000000000
2215#define subst @[y.G@]
2216
2217@ For example, we can use the |IND_GRAPH| feature to create a
2218``wheel'' of $n$ vertices arranged cyclically, all connected to one or
2219more center points. In the directed case, the arcs will run from the
2220center(s) to a cycle; in the undirected case, the edges will join the
2221center(s) to a circuit.
2222
2223@<Applications...@>=
2224Graph *wheel(n,n1,directed)
2225  unsigned long n; /* size of the rim */
2226  unsigned long n1; /* number of center points */
2227  long directed; /* should all arcs go from center to rim and around? */
2228{@+Graph *new_graph=board(2L,0L,0L,0L,1L,0L,directed);
2229                                         /* trivial 2-vertex graph */
2230  if (new_graph) {
2231    new_graph->vertices->ind=n1;
2232    (new_graph->vertices+1)->ind=IND_GRAPH;
2233    (new_graph->vertices+1)->subst=board(n,0L,0L,0L,1L,1L,directed);
2234                             /* cycle or circuit */
2235    new_graph=induced(new_graph,NULL,0L,0L,directed);
2236    if (new_graph) {
2237      sprintf(new_graph->id,"wheel(%lu,%lu,%d)",@|
2238                              n,n1,directed?1:0);
2239    }
2240  }
2241  return new_graph;
2242}
2243
2244@ @(gb_basic.h@>=
2245extern Graph *bi_complete();
2246extern Graph *wheel(); /* standard applications of |induced| */
2247
2248@ @<Basic subroutines@>=
2249Graph *induced(g,description,self,multi,directed)
2250  Graph *g; /* graph marked for induction in its |ind| fields */
2251  char *description; /* string to be mentioned in |new_graph->id| */
2252  long self; /* should self-loops be permitted? */
2253  long multi; /* should multiple arcs be permitted? */
2254  long directed; /* should the graph be directed? */
2255{@+@<Vanilla local variables@>@;@+@q{@>@t}\6{@>@q}@>
2256  register Vertex *u;
2257  register long n=0; /* total number of vertices in induced graph */
2258  register long nn=0; /* number of negative vertices in induced graph */
2259  if (g==NULL) panic(missing_operand); /* where is |g|? */
2260  @<Set up a graph with the induced vertices@>;
2261  @<Insert arcs or edges for induced vertices@>;
2262  @<Restore |g| to its original state@>;
2263  if (gb_trouble_code) {
2264    gb_recycle(new_graph);
2265    panic(alloc_fault); /* aargh, we ran out of memory somewhere back there */
2266  }
2267  return new_graph;
2268}
2269
2270@ @<Set up a graph with the induced vertices@>=
2271@<Determine |n| and |nn|@>;
2272new_graph=gb_new_graph(n);
2273if (new_graph==NULL)
2274  panic(no_room); /* out of memory before we're even started */
2275@<Assign names to the new vertices, and create a map from |g| to |new_graph|@>;
2276sprintf(buffer,",%s,%d,%d,%d)",@|description?description:null_string,@|
2277                    self?1:0,multi?1:0,directed?1:0);
2278make_compound_id(new_graph,"induced(",g,buffer);
2279
2280@ @<Determine |n| and |nn|@>=
2281for (v=g->vertices;v<g->vertices+g->n;v++)
2282  if (v->ind>0) {
2283    if (n>IND_GRAPH) panic(very_bad_specs); /* way too big */
2284    if (v->ind>=IND_GRAPH) {
2285      if (v->subst==NULL) panic(missing_operand+1);
2286        /* substitute graph is missing */
2287      n+=v->subst->n;
2288    }@+else n+=v->ind;
2289  }@+else if (v->ind<-nn) nn=-(v->ind);
2290if (n>IND_GRAPH || nn>IND_GRAPH) panic(very_bad_specs+1); /* gigantic */
2291n+=nn;
2292
2293@ The negative vertices get the negative number as their name. Split vertices
2294get names with an optional prime appended, if the |ind| field is 2;
2295otherwise split vertex names are obtained by appending a colon and an index
2296number between |0| and~|ind-1|. The name of a vertex within a
2297graph |v->subst| is composed of the name of |v| followed by a
2298colon and the name within that graph.
2299
2300We store the original |ind| field in the |mult| field of the first
2301corresponding vertex in the new graph, and change |ind| to point to
2302that vertex. This convention makes it easy
2303to determine the location of each vertex's clone or clones.
2304Of course, if the original |ind| field is zero, we leave it zero (|NULL|),
2305because it has no corresponding vertex in the new graph.
2306
2307@<Assign names to the new vertices, and create a map from |g| to |new_graph|@>=
2308for (k=1,u=new_graph->vertices;k<=nn;k++,u++) {
2309  u->mult=-k;
2310  sprintf(buffer,"%ld",-k);
2311  u->name=gb_save_string(buffer);
2312}
2313for (v=g->vertices;v<g->vertices+g->n;v++)
2314  if ((k=v->ind)<0) v->map=(new_graph->vertices)-(k+1);
2315  else if (k>0) {
2316    u->mult=k;
2317    v->map=u;
2318    if (k<=2) {
2319      u->name=gb_save_string(v->name);
2320      u++;
2321      if (k==2) {
2322        sprintf(buffer,"%s'",v->name);
2323        u->name=gb_save_string(buffer);
2324        u++;
2325      }
2326    }@+else if (k>=IND_GRAPH) @<Make names and arcs for a substituted graph@>@;
2327    else for (j=0;j<k;j++,u++) {
2328      sprintf(buffer,"%.*s:%ld",BUF_SIZE-12,v->name,j);
2329      u->name=gb_save_string(buffer);
2330    }
2331  }
2332
2333@ @<Restore |g| to its original state@>=
2334for (v=g->vertices;v<g->vertices+g->n;v++)
2335  if (v->map) v->ind=v->map->mult;
2336for (v=new_graph->vertices;v<new_graph->vertices+n;v++)
2337  v->u.I=v->v.I=v->z.I=0; /* clear |tmp|, |mult|, |tlen| */
2338
2339@ The heart of the procedure to construct an induced graph is, of course,
2340the part where we map the arcs of |g| into arcs of |new_graph|.
2341
2342Notice that if |v| has a self-loop
2343in the original graph and if |v| is being split into several vertices,
2344it will produce arcs between different clones of itself, but it will not
2345produce self-loops unless |self!=0|. In an undirected graph, a loop
2346from a vertex to itself will not produce multiple edges among its clones,
2347even if |multi!=0|.
2348
2349More precisely, if |v| has |k| clones |u| through |u+k-1|, an original
2350directed arc from |v| to~|v| will generate all $k^2$ possible arcs between
2351them, except that the |k| self-loops will be eliminated when
2352|self==0|.  An original undirected edge from |v| to~|v| will generate
2353$k\choose2$ edges between distinct clones, together with |k|
2354undirected self-loops if |self!=0|.
2355
2356@<Insert arcs or edges for induced vertices@>=
2357for (v=g->vertices;v<g->vertices+g->n;v++) {
2358  u=v->map;
2359  if (u) {@+register Arc *a;@+register Vertex *uu,*vv;
2360    k=u->mult;
2361    if (k<0) k=1; /* |k| is the number of clones of |v| */
2362    else if (k>=IND_GRAPH) k=v->subst->n;
2363    for (;k;k--,u++) {
2364      if (!multi)
2365        @<Take note of existing edges that touch |u|@>;
2366      for (a=v->arcs;a;a=a->next) {
2367        vv=a->tip;
2368        uu=vv->map;
2369        if (uu==NULL) continue;
2370        j=uu->mult;
2371        if (j<0) j=1; /* |j| is the number of clones of |vv| */
2372        else if (j>=IND_GRAPH) j=vv->subst->n;
2373        if (!directed) {
2374          if (vv<v) continue;
2375          if (vv==v) {
2376            if (a->next==a+1) a++; /* skip second half of self-loop */
2377            j=k,uu=u; /* also skip duplicate edges generated by self-loop */
2378          }
2379        }
2380        @<Insert arcs or edges from vertex |u| to vertices
2381              |uu| through |uu+j-1|@>;
2382      }
2383    }
2384  }
2385}
2386
2387@ Again we use the |tmp| and |tlen| trick of |gunion| to handle
2388multiple arcs. (This trick explains why the code in the previous
2389section tries to generate as many arcs as possible from a single
2390vertex |u|, before changing~|u|.)
2391
2392@<Take note of existing edges that touch |u|@>=
2393for (a=u->arcs;a;a=a->next) {
2394  a->tip->tmp=u;
2395  if (directed || a->tip>u || a->next==a+1) a->tip->tlen=a;
2396  else a->tip->tlen=a+1;
2397}
2398
2399@ @<Insert arcs or edges from vertex |u| to vertices |uu|...@>=
2400for (;j;j--,uu++) {
2401  if (u==uu && !self) continue;
2402  if (uu->tmp==u && !multi)
2403    @<Update the minimum arc length from |u| to |uu|, then |continue|@>;
2404  if (directed) gb_new_arc(u,uu,a->len);
2405  else gb_new_edge(u,uu,a->len);
2406  uu->tmp=u;
2407  uu->tlen=((directed || u<=uu)? u->arcs: uu->arcs);
2408}
2409
2410@ @<Update the minimum arc length from |u| to |uu|, then |continue|@>=
2411{@+register Arc *b=uu->tlen; /* existing arc or edge from |u| to |uu| */
2412  if (a->len<b->len) {
2413    b->len=a->len; /* remember the minimum length */
2414    if (!directed) (b+1)->len=a->len;
2415  }
2416  continue;
2417}
2418
2419@ We have now accumulated enough experience to finish off the one
2420remaining piece of program with ease.
2421
2422@<Make names and arcs for a sub...@>=
2423{@+register Graph *gg=v->subst;
2424  register Vertex *vv=gg->vertices;
2425  register Arc *a;
2426  siz_t delta=((siz_t)u)-((siz_t)vv);
2427  for (j=0;j<v->subst->n;j++,u++,vv++) {
2428    sprintf(buffer,"%.*s:%.*s",BUF_SIZE/2-1,v->name,(BUF_SIZE-1)/2,vv->name);
2429    u->name=gb_save_string(buffer);
2430    for (a=vv->arcs;a;a=a->next) {@+register Vertex *vvv=a->tip;
2431      Vertex *uu=vert_offset(vvv,delta);
2432      if (vvv==vv && !self) continue;
2433      if (uu->tmp==u && !multi) @<Update the minimum arc length...@>;
2434      if (!directed) {
2435        if (vvv<vv) continue;
2436        if (vvv==vv && a->next==a+1) a++; /* skip second half of self-loop */
2437        gb_new_edge(u,uu,a->len);
2438      }@+else gb_new_arc(u,uu,a->len);
2439      uu->tmp=u;
2440      uu->tlen=((directed || u<=uu)? u->arcs: uu->arcs);
2441    }
2442  }
2443}
2444
2445@* Index. As usual, we close with an index that
2446shows where the identifiers of \\{gb\_basic} are defined and used.
2447