Path: blob/master/sage/combinat/designs/block_design.py
4079 views
"""1Block designs.23A module to help with constructions and computations of block4designs and other incidence structures.56A block design is an incidence structure consisting of a set of7points P and a set of blocks B, where each block is considered as a8subset of P. More precisely, a *block design* B is a class of9k-element subsets of P such that the number r of blocks that10contain any point x in P is independent of x, and the number lambda11of blocks that contain any given t-element subset T is independent12of the choice of T (see [1] for more). Such a block design is also13called a t-(v,k,lambda)-design, and v (the number of points), b14(the number of blocks), k, r, and lambda are the parameters of the15design. (In Python, lambda is reserved, so we sometimes use lmbda16or L instead.)1718In Sage, sets are replaced by (ordered) lists and the standard19representation of a block design uses P = [0,1,..., v-1], so a20block design is specified by (v,B).2122This software is released under the terms of the GNU General Public23License, version 2 or above (your choice). For details on24licensing, see the accompanying documentation.2526REFERENCES:2728- [1] Block design from wikipedia,29http://en.wikipedia.org/wiki/Block_design3031- [2] What is a block design?,32http://designtheory.org/library/extrep/html/node4.html (in 'The33External Representation of Block Designs' by Peter J. Cameron, Peter34Dobcsanyi, John P. Morgan, Leonard H. Soicher)3536This is a significantly modified form of the module37block_design.py (version 0.6) written by Peter Dobcsanyi38[email protected]. Thanks go to Robert Miller for lots of good39design suggestions.4041Copyright 2007-2008 by Peter Dobcsanyi [email protected], and42David Joyner [email protected].4344TODO: Implement DerivedDesign, ComplementaryDesign,45Hadamard3Design46"""4748from sage.modules.free_module import VectorSpace49from sage.rings.integer_ring import ZZ50from sage.rings.arith import binomial, integer_floor51from sage.combinat.designs.incidence_structures import IncidenceStructure, IncidenceStructureFromMatrix52from sage.misc.decorators import rename_keyword5354### utility functions -------------------------------------------------------5556def tdesign_params(t, v, k, L):57"""58Return the design's parameters: (t, v, b, r , k, L). Note t must be59given.6061EXAMPLES::6263sage: BD = BlockDesign(7,[[0,1,2],[0,3,4],[0,5,6],[1,3,5],[1,4,6],[2,3,6],[2,4,5]])64sage: from sage.combinat.designs.block_design import tdesign_params65sage: tdesign_params(2,7,3,1)66(2, 7, 7, 3, 3, 1)67"""68x = binomial(v, t)69y = binomial(k, t)70b = divmod(L * x, y)[0]71x = binomial(v-1, t-1)72y = binomial(k-1, t-1)73r = integer_floor(L * x/y)74return (t, v, b, r, k, L)7576@rename_keyword(deprecated='Sage version 4.6', method="algorithm")77def ProjectiveGeometryDesign(n, d, F, algorithm=None):78"""79Input: n is the projective dimension, so the number of points is v80= PPn(GF(q)) d is the dimension of the subspaces of P = PPn(GF(q))81which make up the blocks, so b is the number of d-dimensional82subspaces of P8384Wraps GAP Design's PGPointFlatBlockDesign. Does *not* require85GAP's Design.8687EXAMPLES::8889sage: ProjectiveGeometryDesign(2, 1, GF(2))90Incidence structure with 7 points and 7 blocks91sage: BD = ProjectiveGeometryDesign(2, 1, GF(2), algorithm="gap") # requires optional gap package 'design'92sage: BD.is_block_design() # requires optional gap package 'design'93(True, [2, 7, 3, 1])94"""95q = F.order()96from sage.interfaces.gap import gap, GapElement97from sage.sets.set import Set98if algorithm == None:99V = VectorSpace(F, n+1)100points = list(V.subspaces(1))101flats = list(V.subspaces(d+1))102blcks = []103for p in points:104b = []105for i in range(len(flats)):106if p.is_subspace(flats[i]):107b.append(i)108blcks.append(b)109v = (q**(n+1)-1)/(q-1)110return BlockDesign(v, blcks, name="ProjectiveGeometryDesign")111if algorithm == "gap": # Requires GAP's Design112gap.eval('LoadPackage("design")')113gap.eval("D := PGPointFlatBlockDesign( %s, %s, %d )"%(n,q,d))114v = eval(gap.eval("D.v"))115gblcks = eval(gap.eval("D.blocks"))116gB = []117for b in gblcks:118gB.append([x-1 for x in b])119return BlockDesign(v, gB, name="ProjectiveGeometryDesign")120121def AffineGeometryDesign(n, d, F):122r"""123Input: n is the Euclidean dimension, so the number of points is124`v = |F^n|` (F = GF(q), some q) d is the dimension of the125(affine) subspaces of `P = GF(q)^n` which make up the126blocks.127128`AG_{n,d} (F)`, as it is sometimes denoted, is a129`2` - `(v, k, \lambda)` design of points and130`d`- flats (cosets of dimension n) in the affine geometry131`AG_n (F)`, where132133.. math::134135v = q^n,\ k = q^d ,136\lambda =\frac{(q^{n-1}-1) \cdots (q^{n+1-d}-1)}{(q^{n-1}-1) \cdots (q-1)}.137138139140Wraps some functions used in GAP Design's PGPointFlatBlockDesign.141Does *not* require GAP's Design.142143EXAMPLES::144145sage: BD = AffineGeometryDesign(3, 1, GF(2))146sage: BD.parameters()147(2, 8, 2, 2)148sage: BD.is_block_design()149(True, [2, 8, 2, 2])150sage: BD = AffineGeometryDesign(3, 2, GF(2))151sage: BD.parameters()152(2, 8, 4, 12)153sage: BD.is_block_design()154(True, [3, 8, 4, 4])155"""156q = F.order()157from sage.interfaces.gap import gap, GapElement158from sage.sets.set import Set159gap.eval("V:=GaloisField(%s)^%s"%(q,n))160gap.eval("points:=AsSet(V)")161gap.eval("Subs:=AsSet(Subspaces(V,%s));"%d)162gap.eval("CP:=Cartesian(points,Subs)")163flats = gap.eval("flats:=List(CP,x->Sum(x))") # affine spaces164gblcks = eval(gap.eval("AsSortedList(List(flats,f->Filtered([1..Length(points)],i->points[i] in f)));"))165v = q**n166gB = []167for b in gblcks:168gB.append([x-1 for x in b])169return BlockDesign(v, gB, name="AffineGeometryDesign")170171def WittDesign(n):172"""173Input: n is in 9,10,11,12,21,22,23,24.174175Wraps GAP Design's WittDesign. If n=24 then this function returns176the large Witt design W24, the unique (up to isomorphism)1775-(24,8,1) design. If n=12 then this function returns the small178Witt design W12, the unique (up to isomorphism) 5-(12,6,1) design.179The other values of n return a block design derived from these.180181REQUIRES: GAP's Design package.182183EXAMPLES::184185sage: BD = WittDesign(9) # requires optional gap package 'design'186sage: BD.parameters() # requires optional gap package 'design'187(2, 9, 3, 1)188sage: BD # requires optional gap package 'design'189Incidence structure with 9 points and 12 blocks190sage: print BD # requires optional gap package 'design'191WittDesign<points=[0, 1, 2, 3, 4, 5, 6, 7, 8], blocks=[[0, 1, 7], [0, 2, 5], [0, 3, 4], [0, 6, 8], [1, 2, 6], [1, 3, 5], [1, 4, 8], [2, 3, 8], [2, 4, 7], [3, 6, 7], [4, 5, 6], [5, 7, 8]]>192sage: BD = WittDesign(12) # requires optional gap package 'design'193sage: BD.parameters(t=5) # requires optional gap package 'design'194(5, 12, 6, 1)195"""196from sage.interfaces.gap import gap, GapElement197gap.eval('LoadPackage("design")')198gap.eval("B:=WittDesign(%s)"%n)199v = eval(gap.eval("B.v"))200gblcks = eval(gap.eval("B.blocks"))201gB = []202for b in gblcks:203gB.append([x-1 for x in b])204return BlockDesign(v, gB, name="WittDesign", test=True)205206def HadamardDesign(n):207"""208As described in Section 1, p. 10, in [CvL]. The input n must have the209property that there is a Hadamard matrix of order n+1 (and that a210construction of that Hadamard matrix has been implemented...).211212EXAMPLES::213214sage: HadamardDesign(7)215Incidence structure with 7 points and 7 blocks216sage: print HadamardDesign(7)217HadamardDesign<points=[0, 1, 2, 3, 4, 5, 6], blocks=[[0, 1, 2], [0, 3, 4], [0, 5, 6], [1, 3, 5], [1, 4, 6], [2, 3, 6], [2, 4, 5]]>218219REFERENCES:220221- [CvL] P. Cameron, J. H. van Lint, Designs, graphs, codes and222their links, London Math. Soc., 1991.223"""224from sage.combinat.matrices.hadamard_matrix import hadamard_matrix225from sage.matrix.constructor import matrix226H = hadamard_matrix(n+1)227H1 = H.matrix_from_columns(range(1,n+1))228H2 = H1.matrix_from_rows(range(1,n+1))229J = matrix(ZZ,n,n,[1]*n*n)230MS = J.parent()231A = MS((H2+J)/2) # convert -1's to 0's; coerce entries to ZZ232# A is the incidence matrix of the block design233return IncidenceStructureFromMatrix(A,name="HadamardDesign")234235def steiner_triple_system(n):236r"""237Returns a Steiner Triple System238239A Steiner Triple System (STS) of a set `\{0,...,n-1\}`240is a family `S` of 3-sets such that for any `i \not = j`241there exists exactly one set of `S` in which they are242both contained.243244It can alternatively be thought of as a factorization of245the complete graph `K_n` with triangles.246247A Steiner Triple System of a `n`-set exists if and only if248`n \equiv 1 mod 6` or `n \equiv 3 mod 6`, in which case249one can be found through Bose's and Skolem's constructions,250respectively [AndHon97]_.251252INPUT:253254- ``n`` returns a Steiner Triple System of `\{0,...,n-1\}`255256EXAMPLE:257258A Steiner Triple System on `9` elements ::259260sage: from sage.combinat.designs.block_design import steiner_triple_system261sage: sts = steiner_triple_system(9)262sage: sts263Incidence structure with 9 points and 12 blocks264sage: list(sts)265[[0, 1, 5], [0, 2, 4], [0, 3, 6], [0, 7, 8], [1, 2, 3], [1, 4, 7], [1, 6, 8], [2, 5, 8], [2, 6, 7], [3, 4, 8], [3, 5, 7], [4, 5, 6]]266267As any pair of vertices is covered once, its parameters are ::268269sage: sts.parameters()270(2, 9, 3, 1)271272An exception is raised for invalid values of ``n`` ::273274sage: steiner_triple_system(10)275Traceback (most recent call last):276...277ValueError: Steiner triple systems only exist for n = 1 mod 6 or n = 3 mod 6278279REFERENCE:280281.. [AndHon97] A short course in Combinatorial Designs,282Ian Anderson, Iiro Honkala,283Internet Editions, Spring 1997,284http://www.utu.fi/~honkala/designs.ps285"""286287name = "Steiner Triple System on "+str(n)+" elements"288289if n%6 == 3:290t = (n-3)/6291Z = range(2*t+1)292293T = lambda (x,y) : x + (2*t+1)*y294295sts = [[(i,0),(i,1),(i,2)] for i in Z] + \296[[(i,k),(j,k),(((t+1)*(i+j)) % (2*t+1),(k+1)%3)] for k in range(3) for i in Z for j in Z if i != j]297298elif n%6 == 1:299300t = (n-1)/6301N = range(2*t)302T = lambda (x,y) : x+y*t*2 if (x,y) != (-1,-1) else n-1303304L1 = lambda i,j : (i+j) % (int((n-1)/3))305L = lambda i,j : L1(i,j)/2 if L1(i,j)%2 == 0 else t+(L1(i,j)-1)/2306307sts = [[(i,0),(i,1),(i,2)] for i in range(t)] + \308[[(-1,-1),(i,k),(i-t,(k+1) % 3)] for i in range(t,2*t) for k in [0,1,2]] + \309[[(i,k),(j,k),(L(i,j),(k+1) % 3)] for k in [0,1,2] for i in N for j in N if i < j]310311else:312raise ValueError("Steiner triple systems only exist for n = 1 mod 6 or n = 3 mod 6")313314from sage.sets.set import Set315sts = Set(map(lambda x: Set(map(T,x)),sts))316317return BlockDesign(n, sts, name=name)318319def BlockDesign(max_pt, blks, name=None, test=True):320"""321Returns an instance of the IncidenceStructure class. Requires each322B in blks to be contained in range(max_pt). Does not test if the323result is a block design.324325EXAMPLES::326327sage: BlockDesign(7,[[0,1,2],[0,3,4],[0,5,6],[1,3,5],[1,4,6],[2,3,6],[2,4,5]], name="Fano plane")328Incidence structure with 7 points and 7 blocks329sage: print BlockDesign(7,[[0,1,2],[0,3,4],[0,5,6],[1,3,5],[1,4,6],[2,3,6],[2,4,5]], name="Fano plane")330Fano plane<points=[0, 1, 2, 3, 4, 5, 6], blocks=[[0, 1, 2], [0, 3, 4], [0, 5, 6], [1, 3, 5], [1, 4, 6], [2, 3, 6], [2, 4, 5]]>331"""332nm = name333tst = test334if nm == None and test:335nm = "BlockDesign"336BD = BlockDesign_generic( range(max_pt), blks, name=nm, test=tst )337if not(test):338return BD339else:340pars = BD.parameters()341if BD.block_design_checker(pars[0],pars[1],pars[2],pars[3]):342return BD343else:344raise TypeError("parameters are not those of a block design.")345346BlockDesign_generic = IncidenceStructure347"""348Possibly in the future there will be methods which apply to349block designs and not incidence structures. None have been350implemented yet though. The class name BlockDesign_generic351is reserved in the name space in case more specialized352methods are implemented later. In that case, BlockDesign_generic353should inherit from IncidenceStructure.354"""355356357