This repository has been archived by the owner on Jan 30, 2023. It is now read-only.
-
-
Notifications
You must be signed in to change notification settings - Fork 7
/
database.py
3925 lines (3030 loc) · 148 KB
/
database.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
r"""
Database of small combinatorial designs
This module implements known constructions combinatorial designs, and in
particular of :mod:`Orthogonal Arrays
<sage.combinat.designs.orthogonal_arrays>`. Most of them come from the chapter
on :mod:`Mutually Orthogonal Latin Squares
<sage.combinat.designs.latin_squares>` from the Handbook of Combinatorial
Designs.
Most of this would only be a dream without the mathematical knowledge and help
of Julian R. Abel.
All the designs returned by these functions can be obtained through the
``designs.<tab>`` functions.
Implemented constructions :
- :func:`OA(7,18) <OA_7_18>`,
:func:`OA(6,20) <OA_6_20>`,
:func:`OA(7,21) <OA_7_21>`,
:func:`OA(5,22) <OA_5_22>`,
:func:`OA(9,24) <OA_9_24>`,
:func:`OA(6,26) <OA_6_26>`,
:func:`OA(7,28) <OA_7_28>`,
:func:`OA(6,30) <OA_6_30>`,
:func:`OA(7,33) <OA_7_33>`,
:func:`OA(6,34) <OA_6_34>`,
:func:`OA(7,35) <OA_7_35>`,
:func:`OA(10,36) <OA_10_36>`,
:func:`OA(6,38) <OA_6_38>`,
:func:`OA(7,39) <OA_7_39>`,
:func:`OA(9,40) <OA_9_40>`,
:func:`OA(7,42) <OA_7_42>`,
:func:`OA(7,44) <OA_7_44>`,
:func:`OA(8,45) <OA_8_45>`,
:func:`OA(6,46) <OA_6_46>`,
:func:`OA(10,48) <OA_10_48>`,
:func:`OA(8,50) <OA_8_50>`,
:func:`OA(7,51) <OA_7_51>`,
:func:`OA(7,52) <OA_7_52>`,
:func:`OA(7,54) <OA_7_54>`,
:func:`OA(8,55) <OA_8_55>`,
:func:`OA(9,56) <OA_9_56>`,
:func:`OA(9,57) <OA_9_57>`,
:func:`OA(7,60) <OA_7_60>`,
:func:`OA(7,62) <OA_7_62>`,
:func:`OA(7,66) <OA_7_66>`,
:func:`OA(7,68) <OA_7_68>`,
:func:`OA(8,69) <OA_8_69>`,
:func:`OA(7,74) <OA_7_74>`,
:func:`OA(9,75) <OA_9_75>`,
:func:`OA(8,76) <OA_8_76>`,
:func:`OA(11,80) <OA_11_80>`,
:func:`OA(10,82) <OA_10_82>`,
:func:`OA(10,100) <OA_10_100>`,
:func:`OA(15,112) <OA_15_112>`,
:func:`OA(9,120) <OA_9_120>`,
:func:`OA(9,135) <OA_9_135>`,
:func:`OA(12,144) <OA_12_144>`,
:func:`OA(10,154) <OA_10_154>`,
:func:`OA(11,160) <OA_11_160>`,
:func:`OA(16,176) <OA_16_176>`,
:func:`OA(16,208) <OA_16_208>`,
:func:`OA(12,210) <OA_12_210>`,
:func:`OA(15,224) <OA_15_224>`,
:func:`OA(10,262) <OA_10_262>`,
:func:`OA(18,273) <OA_18_273>`,
:func:`OA(12,276) <OA_12_276>`,
:func:`OA(12,298) <OA_12_298>`,
:func:`OA(12,342) <OA_12_342>`,
:func:`OA(20,352) <OA_20_352>`,
:func:`OA(20,416) <OA_20_416>`,
:func:`OA(12,474) <OA_12_474>`,
:func:`OA(9,514) <OA_9_514>`,
:func:`OA(20,544) <OA_20_544>`,
:func:`OA(11,640) <OA_11_640>`,
:func:`OA(10,796) <OA_10_796>`,
:func:`OA(15,896) <OA_15_896>`,
:func:`OA(14,950) <OA_14_950>`,
:func:`OA(33,993) <OA_33_993>`,
:func:`OA(10,1046) <OA_10_1046>`,
:func:`OA(10,1059) <OA_10_1059>`,
:func:`OA(11,2164) <OA_11_2164>`,
:func:`OA(12,3992) <OA_12_3992>`,
:func:`OA(12,3994) <OA_12_3994>`
- :func:`two MOLS of order 10 <MOLS_10_2>`,
:func:`five MOLS of order 12 <MOLS_12_5>`,
:func:`four MOLS of order 14 <MOLS_14_4>`,
:func:`four MOLS of order 15 <MOLS_15_4>`,
:func:`three MOLS of order 18 <MOLS_18_3>`
- :func:`DF(21,5,1) <CDF_21_5_1>`
:func:`DF(25,4,1) <ADF_5x5_4_1>`
:func:`DF(37,4,1) <CDF_37_4_1>`
:func:`DF(81,5,1) <CDF_81_5_1>`
:func:`DF(91,6,1) <CDF_91_6_1>`
:func:`DF(121,5,1) <CDF_121_5_1>`
:func:`DF(141,5,1) <CDF_141_5_1>`
:func:`DF(161,5,1) <CDF_161_5_1>`
:func:`DF(201,5,1) <CDF_201_5_1>`
:func:`DF(221,5,1) <CDF_221_5_1>`
- :func:`RBIBD(120,8,1) <RBIBD_120_8_1>`
**Dictionaries**
The functions defined here are used by
:func:`~sage.combinat.designs.orthogonal_arrays.orthogonal_array`. Thus, the
functions are indexed by dictionary which associates to every integer ``n`` a
pair ``(k,f)`` where ``f`` is a function such that ``f()`` is a `OA(k,n)`. This
dictionary is defined right after the constructions of OA in the file.
The same goes for the constructions of MOLS, used by
:func:`~sage.combinat.designs.latin_squares.mutually_orthogonal_latin_squares`.
REFERENCES:
.. [DesignHandbook] Handbook of Combinatorial Designs (2ed)
Charles Colbourn, Jeffrey Dinitz
Chapman & Hall/CRC
2012
Functions
---------
"""
from sage.combinat.designs.orthogonal_arrays import (OA_from_quasi_difference_matrix,
OA_from_Vmt,
QDM_from_Vmt,
OA_from_wider_OA,
OA_from_PBD,
orthogonal_array)
# Cyclic shift of a list
cyclic_shift = lambda l,i : l[-i:]+l[:-i]
def TD_6_12():
r"""
Returns a `TD(6,12)` as built in [Hanani75]_.
This design is Lemma 3.21 from [Hanani75]_.
EXAMPLE::
sage: from sage.combinat.designs.database import TD_6_12
sage: from sage.combinat.designs.orthogonal_arrays import is_transversal_design
sage: TD = TD_6_12()
sage: is_transversal_design(TD,6,12)
True
The design is available from the general constructor::
sage: designs.transversal_design(6,12,existence=True)
True
REFERENCES:
.. [Hanani75] Haim Hanani,
Balanced incomplete block designs and related designs,
http://dx.doi.org/10.1016/0012-365X(75)90040-0,
Discrete Mathematics, Volume 11, Issue 3, 1975, Pages 255-369.
"""
from sage.groups.additive_abelian.additive_abelian_group import AdditiveAbelianGroup
G = AdditiveAbelianGroup([2,6])
d = [[(0,0),(0,0),(0,0),(0,0),(0,0),(0,0)],
[(0,0),(0,1),(1,0),(0,3),(1,2),(0,4)],
[(0,0),(0,2),(1,2),(1,0),(0,1),(1,5)],
[(0,0),(0,3),(0,2),(0,1),(1,5),(1,4)],
[(0,0),(0,4),(1,1),(1,3),(0,5),(0,2)],
[(0,0),(0,5),(0,1),(1,5),(1,3),(1,1)],
[(0,0),(1,0),(1,3),(0,2),(0,3),(1,2)],
[(0,0),(1,1),(1,5),(1,2),(1,4),(1,0)],
[(0,0),(1,2),(0,4),(0,5),(0,2),(1,3)],
[(0,0),(1,3),(1,4),(0,4),(1,1),(0,1)],
[(0,0),(1,4),(0,5),(1,1),(1,0),(0,3)],
[(0,0),(1,5),(0,3),(1,4),(0,4),(0,5)]]
r = lambda x : int(x[0])*6+int(x[1])
TD = [[i*12+r(G(x)+g) for i,x in enumerate(X)] for X in d for g in G]
for x in TD: x.sort()
return TD
def _MOLS_from_string(s,k):
r"""
Returns MOLS from a string
INPUT:
- ``s`` (string) -- represents the MOLS with entries in a-z. To understand
how the string should be formatted, read the source code of a constructor
that uses it.
- ``k`` (integer) -- the number of MOLS encoded by the string.
EXAMPLES::
sage: _ = designs.mutually_orthogonal_latin_squares(2,10) # indirect doctest
"""
from sage.matrix.constructor import Matrix
matrices = [[] for _ in range(k)]
for i,l in enumerate(s.split()):
l = [ord(x) - 97 for x in l]
matrices[i%k].append(l)
return map(Matrix, matrices)
def MOLS_10_2():
r"""
Returns a pair of MOLS of order 10
Data obtained from
`<http://www.cecm.sfu.ca/organics/papers/lam/paper/html/POLS10/POLS10.html>`_
EXAMPLES::
sage: from sage.combinat.designs.latin_squares import are_mutually_orthogonal_latin_squares
sage: from sage.combinat.designs.database import MOLS_10_2
sage: MOLS = MOLS_10_2()
sage: print are_mutually_orthogonal_latin_squares(MOLS)
True
The design is available from the general constructor::
sage: designs.mutually_orthogonal_latin_squares(2,10,existence=True)
True
"""
from sage.matrix.constructor import Matrix
return [Matrix([[1,8,9,0,2,4,6,3,5,7],
[7,2,8,9,0,3,5,4,6,1],
[6,1,3,8,9,0,4,5,7,2],
[5,7,2,4,8,9,0,6,1,3],
[0,6,1,3,5,8,9,7,2,4],
[9,0,7,2,4,6,8,1,3,5],
[8,9,0,1,3,5,7,2,4,6],
[2,3,4,5,6,7,1,8,9,0],
[3,4,5,6,7,1,2,0,8,9],
[4,5,6,7,1,2,3,9,0,8]]),
Matrix([[1,7,6,5,0,9,8,2,3,4],
[8,2,1,7,6,0,9,3,4,5],
[9,8,3,2,1,7,0,4,5,6],
[0,9,8,4,3,2,1,5,6,7],
[2,0,9,8,5,4,3,6,7,1],
[4,3,0,9,8,6,5,7,1,2],
[6,5,4,0,9,8,7,1,2,3],
[3,4,5,6,7,1,2,8,0,9],
[5,6,7,1,2,3,4,0,9,8],
[7,1,2,3,4,5,6,9,8,0]])]
def MOLS_12_5():
r"""
Returns 5 MOLS of order 12
These MOLS have been found by Brendan McKay.
EXAMPLES::
sage: from sage.combinat.designs.latin_squares import are_mutually_orthogonal_latin_squares
sage: from sage.combinat.designs.database import MOLS_12_5
sage: MOLS = MOLS_12_5()
sage: print are_mutually_orthogonal_latin_squares(MOLS)
True
"""
M = """
abcdefghijkl abcdefghijkl abcdefghijkl abcdefghijkl abcdefghijkl
badcfehgjilk ghefklijcdab dcbahgfelkji jilkbadcfehg klijcdabghef
cdabghefklij efghijklabcd lkjidcbahgfe ijklabcdefgh fehgjilkbadc
dcbahgfelkji cdabghefklij ghefklijcdab badcfehgjilk hgfelkjidcba
ijklabcdefgh klijcdabghef efghijklabcd fehgjilkbadc jilkbadcfehg
jilkbadcfehg fehgjilkbadc hgfelkjidcba dcbahgfelkji lkjidcbahgfe
klijcdabghef hgfelkjidcba jilkbadcfehg cdabghefklij dcbahgfelkji
lkjidcbahgfe ijklabcdefgh badcfehgjilk efghijklabcd ghefklijcdab
efghijklabcd jilkbadcfehg fehgjilkbadc lkjidcbahgfe cdabghefklij
fehgjilkbadc dcbahgfelkji cdabghefklij ghefklijcdab badcfehgjilk
ghefklijcdab badcfehgjilk klijcdabghef hgfelkjidcba ijklabcdefgh
hgfelkjidcba lkjidcbahgfe ijklabcdefgh klijcdabghef efghijklabcd
"""
return _MOLS_from_string(M,5)
def MOLS_14_4():
r"""
Returns four MOLS of order 14
These MOLS were shared by Ian Wanless.
EXAMPLES::
sage: from sage.combinat.designs.latin_squares import are_mutually_orthogonal_latin_squares
sage: from sage.combinat.designs.database import MOLS_14_4
sage: MOLS = MOLS_14_4()
sage: print are_mutually_orthogonal_latin_squares(MOLS)
True
The design is available from the general constructor::
sage: designs.mutually_orthogonal_latin_squares(4,14,existence=True)
True
"""
M = """
bjihgkecalnfmd bfmcenidgjhalk bcdefghijklmna bcdefghijklmna
fckjbhledimagn jcgndfalehkbim gnkjdmiclbhaef jflhnkaecmgdib
mgdlkcbafejnih ikdhaegnmfblcj lifhbjemkangcd emkdjbgfnliahc
cnhemldbigfkaj hjlebifkangcmd dalmgnbjehcfik anighmflkbdcej
edabfnmkcjhgli gbkmfcjeliahdn njcaeifhbdgkml kebcajimdgfhln
nfeicgajldkbhm khclngdafmjibe mfbkcdlagnjihe cgnflembihakjd
iagfjdhnkmelcb elbdmahfignkjc aemnhkjdcifblg ilabkdnhfcjegm
dlnkeafimhcjbg ceabkjnihdmgfl hdnikbagmcelfj ljgnihecbamfdk
gemalfihjnbdkc adficlkmjbenhg cgjflhnbiekdam ndmabcjglfeikh
jhfnimgdbkacel liegjdmhnkcfab fkibmagenldhjc mbhiefljadkncg
hkbgajnmeclidf nmjfhkecbaldgi imhlneckdfajgb difjcnkamehgbl
ablchikgnfdmje fankgbljdcimeh klegafdnhjmcbi ghckmlbdeinjaf
licmdbjfhagenk mgialhcbkedjnf jhadicmlfgbekn fajlgidkhncbme
kmjdneclgbihfa dnhjimbgclfeka ebgcjlkfamindh hkemdacngjblfi
"""
return _MOLS_from_string(M,4)
def MOLS_15_4():
r"""
Returns 4 MOLS of order 15.
These MOLS were shared by Ian Wanless.
EXAMPLES::
sage: from sage.combinat.designs.latin_squares import are_mutually_orthogonal_latin_squares
sage: from sage.combinat.designs.database import MOLS_15_4
sage: MOLS = MOLS_15_4()
sage: print are_mutually_orthogonal_latin_squares(MOLS)
True
The design is available from the general constructor::
sage: designs.mutually_orthogonal_latin_squares(4,15,existence=True)
True
"""
M = """
bcdefghijklmnoa bdgiknfcamehjlo bhealiofmdjgcnk blhcmdinejofakg
abcdefghijklmno acehjlogdbnfikm lcifbmjagnekhdo hcmidnejofkagbl
oabcdefghijklmn nbdfikmahecogjl amdjgcnkbhoflie midnjeofkaglbhc
noabcdefghijklm mocegjlnbifdahk fbnekhdolciagmj dnjeokfaglbhmci
mnoabcdefghijkl lnadfhkmocjgebi kgcoflieamdjbhn jeokfalgbhmcind
lmnoabcdefghijk jmobegilnadkhfc olhdagmjfbnekci ekfalgbmhcindjo
klmnoabcdefghij dknacfhjmobelig jamiebhnkgcofld aflgbmhcnidjoek
jklmnoabcdefghi helobdgiknacfmj ekbnjfciolhdagm lbgmhcnidojekaf
ijklmnoabcdefgh kifmacehjlobdgn nflcokgdjamiebh gmchnidojeakflb
hijklmnoabcdefg oljgnbdfikmaceh iogmdalhekbnjfc chndiojeakfblgm
ghijklmnoabcdef iamkhocegjlnbdf djahnebmiflcokg ndioejakfblgcmh
fghijklmnoabcde gjbnliadfhkmoce hekbiofcnjgmdal ioejafkblgcmhdn
efghijklmnoabcd fhkcomjbegilnad miflcjagdokhneb ojafkbglcmhdnie
defghijklmnoabc egildankcfhjmob cnjgmdkbhealiof fakbglchmdnieoj
cdefghijklmnoab cfhjmeboldgikna gdokhnelcifbmja kgblchmdineojfa
"""
return _MOLS_from_string(M,4)
def MOLS_18_3():
r"""
Returns 3 MOLS of order 18.
These MOLS were shared by Ian Wanless.
EXAMPLES::
sage: from sage.combinat.designs.latin_squares import are_mutually_orthogonal_latin_squares
sage: from sage.combinat.designs.database import MOLS_18_3
sage: MOLS = MOLS_18_3()
sage: print are_mutually_orthogonal_latin_squares(MOLS)
True
The design is available from the general constructor::
sage: designs.mutually_orthogonal_latin_squares(3,18,existence=True)
True
"""
M = """
bgejhkmodcnarilpfq beqpodgcflkrjahnim bcdefghijklmnopqra
echfbilnprdokajmqg gcfrqpehdnmlabkioj rbkamfgdehqjinopcl
qfdigcjmohaeplkbnr ehdgarqfibonmkcljp mlbqhifgajdcrenopk
prgejhdbnaikfqmlco jfiehkargqcponldmb hijbcdefgqraklmnop
oqahfbiecpkjlgrnmd hbgjfilkacrdqpomen gderbkamfpclhqjino
dprkigcjfeqlbmhaon kichbgjmlodaerqpnf fgamlbqhiopkjdcren
geqaljhdbofrmcnikp mljdichbngpekfarqo efghijbcdnopqraklm
chfrkmbieqpgandojl onmbejdicphqflgkar amfgderbkinopclhqj
fdigalncjmrqhkoepb dponcfbejaqirgmhlk qhifgamlbrenopkjdc
lnbqcogpakmhdrifej crhapeoqmkbgidnjfl leqjponacbkidfgmhr
kmocrdphqblnieajgf ndaikqfprmlchjeobg kjmcrponhlbqeafgid
rlnpdaeqigcmojfkbh aoekjlrgqhnmdibfpc dqriklponajbcmhfge
jamoqekfrihdnpbglc rkpflbmahdionejcgq nacleqjpomhrbkidfg
abknprflgdjieoqchm ialqgmcnkrejpofbdh onhkjmcrpgidlbqeaf
hkcloqagmnebjfprdi ljkmrhndoiafbqpgce pondqriklfgeajbcmh
nildmprkhjofcbgqae pmblnaioefjkgcrqhd jponacleqdfgmhrbki
iojmenqalfbpgdchrk fqncmokjpegblhdari crponhkjmeafgidlbq
mjpbnforklgcqhedia qgrodnplbjfhcmieka iklpondqrcmhfgeajb
"""
return _MOLS_from_string(M,3)
# Index of the MOLS constructions
#
# Associates to n the pair (k,f) where f() is a function that returns k MOLS of order n
#
# This dictionary is used by designs.mutually_orthogonal_latin_squares(k,n).
MOLS_constructions = {
10 : (2, MOLS_10_2),
12 : (5, MOLS_12_5),
14 : (4, MOLS_14_4),
15 : (4, MOLS_15_4),
18 : (3, MOLS_18_3)
}
def OA_7_18():
r"""
Returns an OA(7,18)
Proved in [JulianAbel13]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_7_18
sage: OA = OA_7_18()
sage: print is_orthogonal_array(OA,7,18,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(7,18,existence=True)
True
"""
M = """
000 100 100 000 100 100 100 000 000 000 100 000
000 020 100 100 000 120 110 110 010 020 010 120
000 100 022 102 112 001 101 120 121 001 020 002
000 002 100 002 102 122 010 111 110 121 021 001
000 021 000 100 020 112 100 021 112 102 102 012
000 000 011 010 100 010 110 122 011 121 120 111
000 100 002 022 011 121 020 122 100 010 112 112
"""
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing as AdditiveCyclic
from sage.categories.cartesian_product import cartesian_product
G = cartesian_product([AdditiveCyclic(2),AdditiveCyclic(3),AdditiveCyclic(3)])
M = [G(map(int,xxx)) for xxx in M.split()]
M = [M[i*12:(i+1)*12] for i in range(7)]
Mb = [[] for _ in range(7)]
for a,b,c,d,e,f,g in zip(*M):
for y in range(3):
Mb[0].append(a + G((0, 0 , 0 )))
Mb[1].append(b + G((0, 0 , y )))
Mb[2].append(c + G((0, y , 0 )))
Mb[3].append(d + G((0, 2*y , y )))
Mb[4].append(e + G((0, 2*y ,2*y)))
Mb[5].append(f + G((0, y ,2*y)))
Mb[6].append(g + G((0, 0 ,2*y)))
M = OA_from_quasi_difference_matrix(Mb,G,add_col=False)
M = M[:len(M)/2] # only develop w.r.t the last two coordinates
return M
def OA_6_20():
r"""
Returns an OA(6,20)
As explained in the Handbook III.3.49 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_6_20
sage: OA = OA_6_20()
sage: print is_orthogonal_array(OA,6,20,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(6,20,existence=True)
True
"""
M=[[None, 7, 13, 1, 16, 9, 2],
[ 0, 1, 15, 7, 17, 6, 14],
[ 0, 11, 10, 11, 5, 4, 3],
[ 7,None, 13, 16, 1, 2, 9],
[ 1, 0, 15, 17, 7, 14, 6],
[ 11, 0, 10, 5, 11, 3, 4]]
Mb=[[],[],[],[],[],[]]
for R in zip(*M):
a,b,c,d,e,f = R
Mb[0].extend([a,b,c])
Mb[1].extend([b,c,a])
Mb[2].extend([c,a,b])
Mb[3].extend([d,f,e])
Mb[4].extend([e,d,f])
Mb[5].extend([f,e,d])
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing as AdditiveCyclic
M = OA_from_quasi_difference_matrix(Mb,AdditiveCyclic(19),add_col=False)
return M
def OA_7_21():
r"""
Returns an OA(7,21)
As explained in the Handbook III.3.50 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_7_21
sage: OA = OA_7_21()
sage: print is_orthogonal_array(OA,7,21,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(7,21,existence=True)
True
"""
M = [[ 8, 17, 20, 2],
[ 9, 16, 4, 15],
[ 11, 5, 10, 6],
[ 14, 1, 3, 13],
[ 18, 19, 12, 7]]
Mb = [[0],[0],[0],[0],[0],[0]]
for a,b,c,d,e in zip(*M):
Mb[0].extend([a,b,c,d,e])
Mb[1].extend([b,c,d,e,a])
Mb[2].extend([c,d,e,a,b])
Mb[3].extend([d,e,a,b,c])
Mb[4].extend([e,a,b,c,d])
Mb[5].extend([0,0,0,0,0])
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing as AdditiveCyclic
M = OA_from_quasi_difference_matrix(Mb,AdditiveCyclic(21))
return M
def OA_5_22():
r"""
Returns an OA(5,22)
As explained in the Handbook III.3.51 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_5_22
sage: OA = OA_5_22()
sage: print is_orthogonal_array(OA,5,22,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(5,22,existence=True)
True
"""
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing as AdditiveCyclic
G = AdditiveCyclic(21)
M = [
[ 1, 13, 18, 3, 16, 19,None],
[ 16, 19, 1, 13, 18, 3, 0],
[ 18, 3, 16, 19, 1, 13, 0],
[ 6, 15, 6, 15, 6, 15, 0],
[ 12, 9, 19, 16, 5, 2, 0],
]
Mb=[[],[],[],[],[]]
for R in zip(*M):
a,b,c,d,e = [G(x) if x is not None else None for x in R]
Mb[0].extend([a,16*c,4*b])
Mb[1].extend([b,None if a is None else 16*a,4*c])
Mb[2].extend([c,16*b,None if a is None else 4*a])
Mb[3].extend([d,16*d+7,4*d+14])
Mb[4].extend([e,16*e+14,4*e+7])
Mb[0].extend([0,0])
Mb[1].extend([7,14])
Mb[2].extend([14,7])
Mb[3].extend([None,0])
Mb[4].extend([0,None])
M = OA_from_quasi_difference_matrix(Mb,G,add_col=False)
return M
def OA_9_24():
r"""
Returns an OA(9,24)
As explained in the Handbook III.3.52 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_9_24
sage: OA = OA_9_24()
sage: print is_orthogonal_array(OA,9,24,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(9,24,existence=True)
True
"""
M = ("0000 0000 0000 0000 0000 0000 0000 0000 0000 0000 0000 0000 "+
"0000 0010 0100 0110 1000 1010 1100 1110 2000 2010 2100 2110 "+
"0000 0011 1001 2110 0111 2011 2111 1000 0100 1100 1101 2010 "+
"0000 1010 1011 2000 1101 2110 0001 0101 2100 2001 0111 1100 "+
"0000 0001 2010 1111 2111 2100 1101 0011 1010 2101 1000 0110 "+
"0000 1000 2001 1011 0100 1100 0110 2101 2111 0010 1111 2011 "+
"0000 1001 0111 2100 2000 0010 1110 2011 1100 1011 0101 2111 "+
"0000 1011 2101 0100 2110 1001 2000 0110 0101 1111 2011 1010 ")
from sage.groups.additive_abelian.additive_abelian_group import AdditiveAbelianGroup
G = AdditiveAbelianGroup([3,2,2,2])
rlabel = {(x%2,x%3):x for x in range(6)}
M = [G([int(c),int(d),rlabel[int(b),int(a)]]) for a,b,c,d in M.split()]
M = [M[i*12:(i+1)*12] for i in range(8)]
Mb = [[] for _ in range(8)]
for a,b,c,d,e,f,g,h in zip(*M):
Mb[0].extend([a, a + G([0,0,rlabel[0,0]])])
Mb[1].extend([b, b + G([0,1,rlabel[0,0]])])
Mb[2].extend([c, c + G([1,0,rlabel[0,0]])])
Mb[3].extend([d, d + G([1,1,rlabel[0,0]])])
Mb[4].extend([e, e + G([0,0,rlabel[1,0]])])
Mb[5].extend([f, f + G([0,1,rlabel[1,0]])])
Mb[6].extend([g, g + G([1,0,rlabel[1,0]])])
Mb[7].extend([h, h + G([1,1,rlabel[1,0]])])
M = OA_from_quasi_difference_matrix(Mb,G)
return M
def OA_6_26():
r"""
Returns an OA(6,26)
As explained in the Handbook III.3.53 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_6_26
sage: OA = OA_6_26()
sage: print is_orthogonal_array(OA,6,26,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(6,26,existence=True)
True
"""
M = [
[None,None,None,None,None],
[ 0, 0, 0, 0, 0],
[ 1, 6, 7, 8, 14],
[ 3, 11, 20, 18, 10],
[ 6, 10, 14, 1, 5],
[ 4, 19, 5, 12, 2],
]
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing as AdditiveCyclic
G = AdditiveCyclic(21)
Mb=[[0],[0],[0],[0],[0],[0]]
for R in zip(*M):
a,b,c,d,e,f = R
Mb[0].extend([a,b,c,d,e,f])
Mb[1].extend([b,c,d,e,f,a])
Mb[2].extend([c,d,e,f,a,b])
Mb[3].extend([d,e,f,a,b,c])
Mb[4].extend([e,f,a,b,c,d])
Mb[5].extend([f,a,b,c,d,e])
M = OA_from_quasi_difference_matrix(Mb,G,add_col = False)
return M
def OA_7_28():
r"""
Returns an OA(7,28)
As explained in the Handbook III.3.54 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_7_28
sage: OA = OA_7_28()
sage: print is_orthogonal_array(OA,7,28,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(7,28,existence=True)
True
"""
z=2
M = [
[(0,0), (z+1,6),(1,1) ,(1,1) ,(1,3) ,(1,4) ,(0,0) ,(1,4), (z,5) ],
[(z,2), (0,0) ,(1,5) ,(z,1) ,(z,2) ,(z,6) ,(z+1,3),(0,0), (z,1) ],
[(z,3), (z+1,4),(0,0) ,(z+1,5),(z+1,2),(z+1,4),(z+1,2),(1,6), (0,0) ],
[(0,5), (z,6) ,(0,5) ,(0,6) ,(z,3) ,(0,0) ,(0,4) ,(1,5), (z+1,4)],
[(0,3), (0,3) ,(z+1,5),(0,0) ,(0,5) ,(z+1,6),(1,1) ,(0,1), (z,3) ],
[(1,3), (0,6) ,(0,6) ,(1,5) ,(0,0) ,(0,3) ,(z+1,6),(z,2), (0,2) ],
]
from sage.groups.additive_abelian.additive_abelian_group import AdditiveAbelianGroup
from sage.modules.free_module_element import free_module_element as vector
G = AdditiveAbelianGroup([2,2,7])
M = [[G(vector([x//2,x%2,y])) for x,y in L] for L in M]
Mb=[[0],[0],[0],[0],[0],[0]]
for R in zip(*M):
a,b,c,d,e,f = R
Mb[0].extend([a,b,c])
Mb[1].extend([b,c,a])
Mb[2].extend([c,a,b])
Mb[3].extend([d,f,e])
Mb[4].extend([e,d,f])
Mb[5].extend([f,e,d])
M = OA_from_quasi_difference_matrix(Mb,G,add_col = True)
return M
def OA_6_30():
r"""
Returns an OA(6,30)
As explained in the Handbook III.3.55 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_6_30
sage: OA = OA_6_30()
sage: print is_orthogonal_array(OA,6,30,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(6,30,existence=True)
True
"""
M = [
[(0,0),None,(0,0),(0,0),(0,0),(0,0),(0,0)],
[(0,0),(0,0),None,(0,4),(0,2),(0,3),(0,1)],
[(0,0),(3,1),(3,0),None,(4,0),(1,0),(2,0)],
[(0,0),(3,0),(0,2),(1,2),None,(0,1),(0,3)],
[(0,0),(3,3),(1,2),(4,2),(2,0),None,(0,4)],
[(0,0),(4,2),(2,4),(0,3),(2,3),(3,2),None]
]
from sage.groups.additive_abelian.additive_abelian_group import AdditiveAbelianGroup
from sage.modules.free_module_element import free_module_element as vector
G = AdditiveAbelianGroup([5,5])
M = [[None if x is None else G(vector(x)) for x in L] for L in M]
Mb=[[],[],[],[],[],[]]
for R in zip(*M):
a,b,c,d,e,f = R
for i in range(5):
Mb[0].append(None if a is None else a+G(vector((i,i))))
Mb[1].append(None if b is None else b+G(vector((2*i,i))))
Mb[2].append(None if c is None else c+G(vector((i,0))))
Mb[3].append(None if d is None else d+G(vector((4*i,0))))
Mb[4].append(None if e is None else e+G(vector((3*i,4*i))))
Mb[5].append(None if f is None else f+G(vector((4*i,4*i))))
M = OA_from_quasi_difference_matrix(Mb,G,add_col = False)
return M
def OA_7_33():
r"""
Returns an OA(7,33)
As explained in the Handbook III.3.56 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_7_33
sage: OA = OA_7_33()
sage: print is_orthogonal_array(OA,7,33,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(7,33,existence=True)
True
"""
M = [
[ 0, 0, 0, 0, 0, 0],
[ 15, 11, 22, 4, 17, 8],
[ 19, 7, 14, 32, 22, 18],
[ 22, 19, 8, 24, 21, 6],
[ 9, 12, 15, 7, 26, 14],
[ 14, 28, 23, 2, 19, 3]
]
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing as AdditiveCyclic
G = AdditiveCyclic(33)
Mb=[[0,1,7],[0,4,28],[0,16,13],[0,31,19],[0,25,10],[0,22,0]]
for R in zip(*M):
a,b,c,d,e,f = R
for i in range(5):
Mb[0].append(a)
Mb[1].append(b)
Mb[2].append(c)
Mb[3].append(d)
Mb[4].append(e)
Mb[5].append(f)
a,b,c,d,e,f = 4*e,4*a,4*b,4*c,4*d,4*f
M = OA_from_quasi_difference_matrix(Mb,G,add_col = True)
return M
def OA_6_34():
r"""
Returns an OA(6,34)
As explained in the Handbook III.3.57 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_6_34
sage: OA = OA_6_34()
sage: print is_orthogonal_array(OA,6,34,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(6,34,existence=True)
True
"""
M = [
[None, 0, 0, 0, 0, 0],
[ 30, 17, 10, 25, 23, 8],
[ 22, 4, 32, 29, 28, 22],
[ 25, 10, 20, 15, 21, 16],
[ 0, 12, 15, 16, 32, 23],
[ 6, 11, 18, 14, 9, 20]
]
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing as AdditiveCyclic
G = AdditiveCyclic(33)
Mb=[[0,1,3,10,5],[0,4,12,7,20],[0,16,15,28,14],[0,31,27,13,23],[0,25,9,19,26],[0,11,11,0,None]]
times4 = lambda x : None if x is None else 4*x
for R in zip(*M):
a,b,c,d,e,f = [None if x is None else G(x) for x in R]
for i in range(5):
Mb[0].append(a)
Mb[1].append(b)
Mb[2].append(c)
Mb[3].append(d)
Mb[4].append(e)
Mb[5].append(f)
a,b,c,d,e,f = map(times4,[e,a,b,c,d,f])
M = OA_from_quasi_difference_matrix(Mb,G,add_col = False)
return M
def OA_7_35():
r"""
Returns an OA(7,35)
As explained in the Handbook III.3.58 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_7_35
sage: OA = OA_7_35()
sage: print is_orthogonal_array(OA,7,35,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(7,35,existence=True)
True
"""
M = [
[ 0, 15, 30, 10, 25, 1, 16, 31, 11, 26, 2, 17, 32, 12, 6, 3, 18, 33, 27, 21, 4, 19, 13, 7, 22, 5, 34, 28, 8, 23, 20, 14, 29, 9, 24],
[ 0, 22, 16, 3, 4, 9, 10, 32, 26, 13, 18, 5, 27, 14, 15, 20, 7, 1, 23, 31, 29, 2, 24, 11, 19, 17, 25, 12, 6, 28, 33, 34, 21, 8, 30],
[ 0, 29, 2, 31, 18, 10, 32, 26, 34, 28, 27, 21, 15, 9, 17, 30, 3, 4, 5, 20, 12, 6, 14, 22, 16, 8, 23, 24, 25, 33, 11, 19, 13, 7, 1],
[ 0, 8, 9, 17, 11, 25, 19, 27, 28, 1, 15, 23, 31, 4, 26, 12, 6, 14, 29, 16, 2, 3, 18, 33, 34, 20, 7, 22, 30, 24, 10, 32, 5, 13, 21],
[ 0, 1, 23, 24, 32, 33, 6, 7, 29, 30, 10, 11, 12, 13, 28, 8, 9, 31, 4, 5, 27, 14, 15, 16, 3, 25, 26, 34, 21, 22, 2, 17, 18, 19, 20],
[0]*35
]
from sage.rings.finite_rings.integer_mod_ring import IntegerModRing as AdditiveCyclic
G = AdditiveCyclic(35)
M = OA_from_quasi_difference_matrix(M,G,add_col = True)
return M
def OA_10_36():
r"""
Returns an OA(10,36)
As explained in the Handbook III.3.59 [DesignHandbook]_.
.. SEEALSO::
:func:`sage.combinat.designs.orthogonal_arrays.OA_from_quasi_difference_matrix`
EXAMPLES::
sage: from sage.combinat.designs.designs_pyx import is_orthogonal_array
sage: from sage.combinat.designs.database import OA_10_36
sage: OA = OA_10_36()
sage: print is_orthogonal_array(OA,10,36,2)
True
The design is available from the general constructor::
sage: designs.orthogonal_array(10,36,existence=True)
True
"""
M = [
[(0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0), (0,0,0,0)],
[(0,0,0,0), (0,1,0,0), (1,0,0,0), (1,1,0,0), (0,0,0,1), (0,1,0,1), (1,0,0,1), (1,1,0,1), (0,0,0,2), (0,1,0,2), (1,0,0,2), (1,1,0,2)],
[(0,0,0,0), (1,1,1,2), (0,0,2,1), (0,0,1,2), (0,1,2,0), (0,1,0,2), (1,1,1,1), (0,1,1,1), (1,1,1,0), (1,0,2,2), (1,0,0,1), (1,0,1,0)],
[(0,0,0,0), (0,0,1,0), (1,0,1,0), (0,1,0,0), (1,1,0,0), (1,0,2,0), (1,0,0,0), (0,1,2,0), (1,1,2,0), (0,0,2,0), (1,1,1,0), (0,1,1,0)],
[(0,0,0,0), (0,1,2,0), (0,0,1,0), (1,1,1,0), (1,0,2,0), (1,0,1,0), (0,1,0,0), (0,0,2,0), (0,1,1,0), (1,1,0,0), (1,1,2,0), (1,0,0,0)],
[(0,0,0,0), (0,1,1,0), (0,1,2,0), (1,1,2,0), (1,1,0,2), (0,0,1,2), (1,1,2,2), (1,0,0,2), (1,0,0,1), (1,0,1,1), (0,0,2,1), (0,1,1,1)],
[(0,0,0,0), (1,0,1,0), (1,1,0,1), (1,0,1,2), (1,0,2,2), (0,0,2,1), (0,1,0,1), (0,1,0,0), (1,1,2,2), (0,1,1,0), (0,0,1,2), (1,1,2,1)],
[(0,0,0,0), (1,1,0,0), (0,1,1,0), (1,0,2,1), (0,1,0,2), (1,0,2,2), (0,0,2,2), (1,1,1,0), (1,0,1,1), (0,1,2,1), (1,1,1,1), (0,0,0,2)],
[(0,0,0,0), (1,0,0,0), (1,1,1,0), (0,1,1,2), (1,1,2,1), (0,1,1,1), (0,0,1,1), (1,0,2,0), (0,1,2,2), (1,1,0,2), (1,0,2,2), (0,0,0,1)]
]
from sage.groups.additive_abelian.additive_abelian_group import AdditiveAbelianGroup
from sage.modules.free_module_element import free_module_element as vector
G = AdditiveAbelianGroup([2,2,3,3])
M = [[G(vector(x)) for x in L] for L in M]
Mb=[[],[],[],[],[],[],[],[],[]]
for R in zip(*M):
a,b,c,d,e,f,g,h,i = R
for y in range(3):
Mb[0].append(a+G(vector([0,0,0,0])))
Mb[1].append(b+G(vector([0,0,y,0])))