Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/combinat/alternating_sign_matrix.py
8817 views
1
r"""
2
Alternating Sign Matrices
3
4
AUTHORS:
5
6
- Mike Hansen (2007): Initial version
7
- Pierre Cange, Luis Serrano (2012): Added monotone triangles
8
- Travis Scrimshaw (2013-28-03): Added element class for ASM's and made
9
:class:`MonotoneTriangles` inherit from :class:`GelfandTsetlinPatterns`
10
- Jessica Striker (2013): Added additional methods
11
"""
12
#*****************************************************************************
13
# Copyright (C) 2007 Mike Hansen <[email protected]>,
14
# 2012 Pierre Cagne <[email protected]>,
15
# Luis Serrano <[email protected]>
16
# 2013 Travis Scrimshaw <[email protected]>
17
# 2013 Jessica Striker <[email protected]>
18
#
19
# Distributed under the terms of the GNU General Public License (GPL)
20
#
21
# This code is distributed in the hope that it will be useful, but
22
# WITHOUT ANY WARRANTY; without even the implied warranty of
23
# MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
24
# General Public License for more details.
25
#
26
# The full text of the GPL is available at:
27
#
28
# http://www.gnu.org/licenses/
29
#*****************************************************************************
30
31
import itertools
32
import copy
33
from sage.misc.classcall_metaclass import ClasscallMetaclass
34
from sage.misc.flatten import flatten
35
from sage.misc.misc import prod
36
from sage.structure.unique_representation import UniqueRepresentation
37
from sage.structure.parent import Parent
38
from sage.structure.element import Element
39
from sage.categories.finite_enumerated_sets import FiniteEnumeratedSets
40
from sage.matrix.matrix_space import MatrixSpace
41
from sage.matrix.constructor import matrix
42
from sage.rings.all import ZZ, factorial
43
from sage.rings.integer import Integer
44
from sage.combinat.posets.lattices import LatticePoset
45
from sage.combinat.gelfand_tsetlin_patterns import GelfandTsetlinPatternsTopRow
46
from sage.sets.set import Set
47
from sage.combinat.combinatorial_map import combinatorial_map
48
from sage.combinat.non_decreasing_parking_function import NonDecreasingParkingFunction
49
from sage.combinat.permutation import Permutation
50
51
class AlternatingSignMatrix(Element):
52
r"""
53
An alternating sign matrix.
54
55
An alternating sign matrix is a square matrix of `0`'s, `1`'s and `-1`'s
56
such that the sum of each row and column is `1` and the non-zero
57
entries in each row and column alternate in sign.
58
"""
59
__metaclass__ = ClasscallMetaclass
60
61
@staticmethod
62
def __classcall_private__(cls, asm):
63
"""
64
Create an ASM.
65
66
EXAMPLES::
67
68
sage: AlternatingSignMatrix([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
69
[1 0 0]
70
[0 1 0]
71
[0 0 1]
72
"""
73
asm = matrix(asm)
74
if not asm.is_square():
75
raise ValueError("The alternating sign matrices must be square")
76
P = AlternatingSignMatrices(asm.nrows())
77
if asm not in P:
78
raise ValueError("Invalid alternating sign matrix")
79
return P(asm)
80
81
def __init__(self, parent, asm):
82
"""
83
Initialize ``self``.
84
85
EXAMPLES:
86
87
sage: A = AlternatingSignMatrices(3)
88
sage: elt = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
89
sage: TestSuite(elt).run()
90
"""
91
self._matrix = asm
92
Element.__init__(self, parent)
93
94
def _repr_(self):
95
"""
96
Return a string representation of ``self``.
97
98
EXAMPLES::
99
100
sage: A = AlternatingSignMatrices(3)
101
sage: A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
102
[1 0 0]
103
[0 1 0]
104
[0 0 1]
105
"""
106
return repr(self._matrix)
107
108
def __eq__(self, other):
109
"""
110
Check equality.
111
112
EXAMPLES::
113
114
sage: A = AlternatingSignMatrices(3)
115
sage: M = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
116
sage: M == A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
117
True
118
sage: M == A([[1, 0, 0],[0, 0, 1],[0, 1, 0]])
119
False
120
"""
121
if isinstance(other, AlternatingSignMatrix):
122
return self._matrix == other._matrix
123
return self._matrix == other
124
125
def __ne__(self, other):
126
"""
127
Check not equals. This is needed, see :trac:`14762`.
128
129
EXAMPLES::
130
131
sage: A = AlternatingSignMatrices(3)
132
sage: M = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
133
sage: M != A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
134
False
135
sage: M != A([[1, 0, 0],[0, 0, 1],[0, 1, 0]])
136
True
137
"""
138
return not self.__eq__(other)
139
140
def __le__(self, other):
141
"""
142
Check less than or equal to. This is needed, see :trac:`15372`.
143
144
EXAMPLES::
145
146
sage: A = AlternatingSignMatrices(3)
147
sage: M = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
148
sage: M <= A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
149
True
150
sage: M <= A([[1, 0, 0],[0, 0, 1],[0, 1, 0]])
151
False
152
"""
153
if isinstance(other, AlternatingSignMatrix):
154
return self._matrix <= other._matrix
155
return False #return False if other is not an ASM
156
157
def __lt__(self, other):
158
"""
159
Check less than. This is needed, see :trac:`15372`.
160
161
EXAMPLES::
162
163
sage: A = AlternatingSignMatrices(3)
164
sage: M = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
165
sage: M < A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
166
False
167
"""
168
if isinstance(other, AlternatingSignMatrix):
169
return self._matrix < other._matrix
170
return False #return False if other is not an ASM
171
172
def __ge__(self, other):
173
"""
174
Check greater than or equal to. This is needed, see :trac:`15372`.
175
176
EXAMPLES::
177
178
sage: A = AlternatingSignMatrices(3)
179
sage: M = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
180
sage: M >= A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
181
True
182
sage: M >= A([[1, 0, 0],[0, 0, 1],[0, 1, 0]])
183
True
184
"""
185
if isinstance(other, AlternatingSignMatrix):
186
return self._matrix >= other._matrix
187
return False #return False if other is not an ASM
188
189
def __gt__(self, other):
190
"""
191
Check greater than. This is needed, see :trac:`15372`.
192
193
EXAMPLES::
194
195
sage: A = AlternatingSignMatrices(3)
196
sage: M = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
197
sage: M > A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
198
False
199
"""
200
if isinstance(other, AlternatingSignMatrix):
201
return self._matrix > other._matrix
202
return False #return False if other is not an ASM
203
204
def _latex_(self):
205
r"""
206
Return a `\LaTeX` representation of ``self``.
207
208
EXAMPLES::
209
210
sage: A = AlternatingSignMatrices(3)
211
sage: latex(A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]))
212
\left(\begin{array}{rrr}
213
1 & 0 & 0 \\
214
0 & 1 & 0 \\
215
0 & 0 & 1
216
\end{array}\right)
217
"""
218
return self._matrix._latex_()
219
220
def to_matrix(self):
221
"""
222
Return ``self`` as a regular matrix.
223
224
EXAMPLES::
225
226
sage: A = AlternatingSignMatrices(3)
227
sage: asm = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]])
228
sage: m = asm.to_matrix(); m
229
[1 0 0]
230
[0 1 0]
231
[0 0 1]
232
sage: m.parent()
233
Full MatrixSpace of 3 by 3 dense matrices over Integer Ring
234
"""
235
return copy.copy(self._matrix)
236
237
@combinatorial_map(name='to monotone triangle')
238
def to_monotone_triangle(self):
239
r"""
240
Return a monotone triangle from ``self``.
241
242
EXAMPLES::
243
244
sage: A = AlternatingSignMatrices(3)
245
sage: A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]).to_monotone_triangle()
246
[[3, 2, 1], [2, 1], [1]]
247
sage: asm = A([[0, 1, 0],[1, -1, 1],[0, 1, 0]])
248
sage: asm.to_monotone_triangle()
249
[[3, 2, 1], [3, 1], [2]]
250
sage: asm = A([[0, 0, 1],[1, 0, 0],[0, 1, 0]])
251
sage: asm.to_monotone_triangle()
252
[[3, 2, 1], [3, 1], [3]]
253
sage: A.from_monotone_triangle(asm.to_monotone_triangle()) == asm
254
True
255
"""
256
n = self._matrix.nrows()
257
triangle = [None]*n
258
prev = [0]*n
259
for j, row in enumerate(self._matrix):
260
add_row = [a+b for (a,b) in itertools.izip(row, prev)]
261
line = [i+1 for (i,val) in enumerate(add_row) if val==1]
262
triangle[n-1-j] = list(reversed(line))
263
prev = add_row
264
return MonotoneTriangles(n)(triangle)
265
266
@combinatorial_map(name='rotate counterclockwise')
267
def rotate_ccw(self):
268
r"""
269
Return the counterclockwise quarter turn rotation of ``self``.
270
271
EXAMPLES::
272
273
sage: A = AlternatingSignMatrices(3)
274
sage: A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]).rotate_ccw()
275
[0 0 1]
276
[0 1 0]
277
[1 0 0]
278
sage: asm = A([[0, 0, 1],[1, 0, 0],[0, 1, 0]])
279
sage: asm.rotate_ccw()
280
[1 0 0]
281
[0 0 1]
282
[0 1 0]
283
"""
284
l = list(self._matrix.transpose())
285
l.reverse()
286
return AlternatingSignMatrix(matrix(l))
287
288
@combinatorial_map(name='rotate clockwise')
289
def rotate_cw(self):
290
r"""
291
Return the clockwise quarter turn rotation of ``self``.
292
293
EXAMPLES::
294
295
sage: A = AlternatingSignMatrices(3)
296
sage: A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]).rotate_cw()
297
[0 0 1]
298
[0 1 0]
299
[1 0 0]
300
sage: asm = A([[0, 0, 1],[1, 0, 0],[0, 1, 0]])
301
sage: asm.rotate_cw()
302
[0 1 0]
303
[1 0 0]
304
[0 0 1]
305
"""
306
l = list(self._matrix.transpose())
307
l.reverse()
308
return AlternatingSignMatrix(matrix(l).transpose().antitranspose())
309
310
@combinatorial_map(name='transpose')
311
def transpose(self):
312
r"""
313
Return the counterclockwise quarter turn rotation of ``self``.
314
315
EXAMPLES::
316
317
sage: A = AlternatingSignMatrices(3)
318
sage: A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]).transpose()
319
[1 0 0]
320
[0 1 0]
321
[0 0 1]
322
sage: asm = A([[0, 0, 1],[1, 0, 0],[0, 1, 0]])
323
sage: asm.transpose()
324
[0 1 0]
325
[0 0 1]
326
[1 0 0]
327
"""
328
return AlternatingSignMatrix(self._matrix.transpose())
329
330
def corner_sum_matrix(self):
331
r"""
332
Return the corner sum matrix from ``self``.
333
334
EXAMPLES::
335
336
sage: A = AlternatingSignMatrices(3)
337
sage: A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]).corner_sum_matrix()
338
[0 0 0 0]
339
[0 1 1 1]
340
[0 1 2 2]
341
[0 1 2 3]
342
sage: asm = A([[0, 1, 0],[1, -1, 1],[0, 1, 0]])
343
sage: asm.corner_sum_matrix()
344
[0 0 0 0]
345
[0 0 1 1]
346
[0 1 1 2]
347
[0 1 2 3]
348
sage: asm = A([[0, 0, 1],[1, 0, 0],[0, 1, 0]])
349
sage: asm.corner_sum_matrix()
350
[0 0 0 0]
351
[0 0 1 1]
352
[0 0 1 2]
353
[0 1 2 3]
354
"""
355
asm = self.to_matrix()
356
n = asm.nrows() + 1
357
return matrix([[nw_corner_sum(asm,i,j) for i in range(n)] for j in range(n)])
358
359
def height_function(self):
360
r"""
361
Return the height function from ``self``. A height function
362
corresponding to an `n \times n` ASM is an `(n+1) \times (n+1)` matrix
363
such that the first row is `0, 1, \ldots, n`, the last row is
364
`n, n-1, \ldots, 1, 0`, and the difference between adjacent entries
365
is 1.
366
367
EXAMPLES::
368
369
sage: A = AlternatingSignMatrices(3)
370
sage: A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]).height_function()
371
[0 1 2 3]
372
[1 0 1 2]
373
[2 1 0 1]
374
[3 2 1 0]
375
sage: asm = A([[0, 1, 0],[1, -1, 1],[0, 1, 0]])
376
sage: asm.height_function()
377
[0 1 2 3]
378
[1 2 1 2]
379
[2 1 2 1]
380
[3 2 1 0]
381
sage: asm = A([[0, 0, 1],[1, 0, 0],[0, 1, 0]])
382
sage: asm.height_function()
383
[0 1 2 3]
384
[1 2 1 2]
385
[2 3 2 1]
386
[3 2 1 0]
387
"""
388
asm = self.to_matrix()
389
n = asm.nrows() + 1
390
return matrix([[i+j-2*nw_corner_sum(asm,i,j) for i in range(n)] for j in range(n)])
391
392
@combinatorial_map(name='gyration')
393
def gyration(self):
394
r"""
395
Return the alternating sign matrix obtained by applying the gyration
396
action to the height function in bijection with ``self``.
397
398
Gyration acts on height functions as follows. Go through the entries of
399
the matrix, first those for which the sum of the row and column indices
400
is even, then for those for which it is odd, and increment or decrement
401
the squares by 2 wherever possible such that the resulting matrix is
402
still a height function. Gyration was first defined in [Wieland00]_ as
403
an action on fully-packed loops.
404
405
REFERENCES:
406
407
.. [Wieland00] B. Wieland. *A large dihedral symmetry of the set of
408
alternating sign matrices*. Electron. J. Combin. 7 (2000).
409
410
EXAMPLES::
411
412
sage: A = AlternatingSignMatrices(3)
413
sage: A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]).gyration()
414
[0 0 1]
415
[0 1 0]
416
[1 0 0]
417
sage: asm = A([[0, 1, 0],[1, -1, 1],[0, 1, 0]])
418
sage: asm.gyration()
419
[1 0 0]
420
[0 1 0]
421
[0 0 1]
422
sage: asm = A([[0, 0, 1],[1, 0, 0],[0, 1, 0]])
423
sage: asm.gyration()
424
[0 1 0]
425
[0 0 1]
426
[1 0 0]
427
"""
428
A = self.parent()
429
hf = list(self.height_function())
430
k = len(hf) - 1
431
for i in range(1,k):
432
for j in range(1,k):
433
if (i+j) % 2 == 0 \
434
and hf[i-1][j] == hf[i+1][j] == hf[i][j+1] == hf[i][j-1]:
435
if hf[i][j] < hf[i+1][j]:
436
hf[i][j] += 2
437
else:
438
hf[i][j] -= 2
439
for i in range(1,k):
440
for j in range(1,k):
441
if (i+j) % 2 == 1 \
442
and hf[i-1][j] == hf[i+1][j] == hf[i][j+1] == hf[i][j-1]:
443
if hf[i][j] < hf[i+1][j]:
444
hf[i][j] += 2
445
else:
446
hf[i][j] -= 2
447
return A.from_height_function(matrix(hf))
448
449
def ASM_compatible(self, B):
450
r"""
451
Return ``True`` if ``self`` and ``B`` are compatible alternating sign
452
matrices in the sense of [EKLP92]_. (If ``self`` is of size `n`, ``B``
453
must be of size `n+1`.)
454
455
In [EKLP92]_, there is a notion of a pair of ASM's with sizes differing
456
by 1 being compatible, in the sense that they can be combined to encode
457
a tiling of the Aztec Diamond.
458
459
REFERENCES:
460
461
.. [EKLP92] N. Elkies, G. Kuperberg, M. Larsen, J. Propp,
462
*Alternating-Sign Matrices and Domino Tilings*, Journal of Algebraic
463
Combinatorics, volume 1 (1992), p. 111-132.
464
465
EXAMPLES::
466
467
sage: A = AlternatingSignMatrix(matrix([[0,0,1,0],[0,1,-1,1],[1,0,0,0],[0,0,1,0]]))
468
sage: B = AlternatingSignMatrix(matrix([[0,0,1,0,0],[0,0,0,1,0],[1,0,0,-1,1],[0,1,0,0,0],[0,0,0,1,0]]))
469
sage: A.ASM_compatible(B)
470
True
471
sage: A = AlternatingSignMatrix(matrix([[0,1,0],[1,-1,1],[0,1,0]]))
472
sage: B = AlternatingSignMatrix(matrix([[0,0,1,0],[0,0,0,1],[1,0,0,0],[0,1,0,0]]))
473
sage: A.ASM_compatible(B)
474
False
475
"""
476
if B.parent()._n - self.parent()._n != 1:
477
raise ValueError("mismatched sizes")
478
479
AA = self.corner_sum_matrix()
480
BB = B.corner_sum_matrix()
481
for i in range(0, len(AA[0])):
482
for j in range(0, len(AA[0])):
483
if not (AA[i,j]>=BB[i,j] and AA[i,j]>=BB[i+1,j+1]-1 \
484
and AA[i,j]<=BB[i+1,j] and AA[i,j]<=BB[i,j+1]):
485
return False
486
return True
487
488
def ASM_compatible_bigger(self):
489
r"""
490
Return all ASM's compatible with ``self`` that are of size one greater
491
than ``self``.
492
493
Given an `n \times n` alternating sign matrix `A`, there are as many
494
ASM's of size `n+1` compatible with `A` as 2 raised to the power of
495
the number of 1's in `A` [EKLP92]_.
496
497
EXAMPLES::
498
499
sage: A = AlternatingSignMatrix(matrix([[1,0],[0,1]]))
500
sage: A.ASM_compatible_bigger()
501
[
502
[ 0 1 0] [1 0 0] [0 1 0] [1 0 0]
503
[ 1 -1 1] [0 0 1] [1 0 0] [0 1 0]
504
[ 0 1 0], [0 1 0], [0 0 1], [0 0 1]
505
]
506
sage: B = AlternatingSignMatrix(matrix([[0,1],[1,0]]))
507
sage: B.ASM_compatible_bigger()
508
[
509
[0 0 1] [0 0 1] [0 1 0] [ 0 1 0]
510
[0 1 0] [1 0 0] [0 0 1] [ 1 -1 1]
511
[1 0 0], [0 1 0], [1 0 0], [ 0 1 0]
512
]
513
"""
514
n = self.parent()._n + 1
515
M = AlternatingSignMatrices(n)
516
sign = []
517
asm = self.to_matrix()
518
B = matrix(n+1)
519
A = matrix([[2*(i+j-2*nw_corner_sum(asm,i,j))+1 for i in range(n)]
520
for j in range(n)])
521
for a in range(n+1):
522
B[a,0] = 2*a
523
B[0,a] = 2*a
524
B[a,n] = 2*(n-a)
525
B[n,a] = 2*(n-a)
526
527
for i in range(1,n):
528
for j in range(1,n):
529
if A[i-1,j-1] == A[i,j] == A[i-1,j]-2 == A[i,j-1]-2:
530
B[i,j] = -A[i,j]
531
sign.append([i,j])
532
else:
533
B[i,j] = list({A[i-1,j-1]-1,A[i-1,j-1]+3} & {A[i-1,j]-3,A[i-1,j]+1} & {A[i,j-1]-3,A[i,j-1]+1} & {A[i,j]-1,A[i,j]+3})[0]
534
535
output = [B]
536
for b in range(len(sign)):
537
N = len(output)
538
for c in range(N):
539
d = copy.copy(output[c])
540
output[c][sign[b][0],sign[b][1]] = -output[c][sign[b][0], sign[b][1]] + 3
541
d[sign[b][0],sign[b][1]] = -d[sign[b][0], sign[b][1]]-1
542
output.append(d)
543
544
for k in range(len(output)):
545
output[k] = M.from_height_function(output[k]/2)
546
return(output)
547
548
def ASM_compatible_smaller(self):
549
r"""
550
Return the list of all ASMs compatible with ``self`` that are of size
551
one smaller than ``self``.
552
553
Given an alternating sign matrix `A` of size `n`, there are as many
554
ASM's of size `n-1` compatible with it as 2 raised to the power of
555
the number of `-1`'s in `A` [EKLP92]_.
556
557
EXAMPLES::
558
559
sage: A = AlternatingSignMatrix(matrix([[0,0,1,0],[0,1,-1,1],[1,0,0,0],[0,0,1,0]]))
560
sage: A.ASM_compatible_smaller()
561
[
562
[0 0 1] [ 0 1 0]
563
[1 0 0] [ 1 -1 1]
564
[0 1 0], [ 0 1 0]
565
]
566
sage: B = AlternatingSignMatrix(matrix([[1,0,0],[0,0,1],[0,1,0]]))
567
sage: B.ASM_compatible_smaller()
568
[
569
[1 0]
570
[0 1]
571
]
572
573
"""
574
n = self.parent()._n
575
M = AlternatingSignMatrices(n)
576
A = matrix(n)
577
asm = self.to_matrix()
578
B = matrix([[2*(i+j-2*nw_corner_sum(asm,i,j)) for i in range(n)] for j in range(n)])
579
sign = []
580
for a in range(n):
581
A[a,0] = 2*a + 1
582
A[0,a] = 2*a + 1
583
A[n-1,a] = 2*(n-a) - 1
584
A[a,n-1] = 2*(n-a) - 1
585
586
for i in range(n-1):
587
for j in range(n-1):
588
if B[i+1,j+1] == B[i,j] == B[i,j+1]+2 == B[i+1,j]+2:
589
A[i,j] = -B[i,j]
590
sign.append([i,j])
591
else:
592
A[i,j] = list({B[i,j]+1,B[i,j]-3} & {B[i,j+1]+3,B[i,j+1]-1} & {B[i+1,j]+3,B[i+1,j]-1} & {B[i+1,j+1]+1,B[i+1,j+1]-3})[0]
593
594
output = [A]
595
for b in range(len(sign)):
596
N = len(output)
597
for c in range(N):
598
d = copy.copy(output[c])
599
output[c][sign[b][0],sign[b][1]] = -output[c][sign[b][0], sign[b][1]]+1
600
d[sign[b][0],sign[b][1]] = -d[sign[b][0], sign[b][1]]-3
601
output.append(d)
602
for k in range(0,len(output)):
603
output[k] = M.from_height_function((output[k]-matrix.ones(n,n))/2)
604
return(output)
605
606
@combinatorial_map(name='to Dyck word')
607
def to_dyck_word(self):
608
r"""
609
Return the Dyck word determined by the last diagonal of
610
the monotone triangle corresponding to ``self``.
611
612
EXAMPLES::
613
614
sage: A = AlternatingSignMatrices(3)
615
sage: A([[0,1,0],[1,0,0],[0,0,1]]).to_dyck_word()
616
[1, 1, 0, 0, 1, 0]
617
sage: d = A([[0,1,0],[1,-1,1],[0,1,0]]).to_dyck_word(); d
618
[1, 1, 0, 1, 0, 0]
619
sage: parent(d)
620
Complete Dyck words
621
"""
622
MT = self.to_monotone_triangle()
623
nplus = self._matrix.nrows() + 1
624
parkfn = [nplus - row[0] for row in list(MT) if len(row) > 0]
625
return NonDecreasingParkingFunction(parkfn).to_dyck_word().reverse()
626
627
def number_negative_ones(self):
628
"""
629
Return the number of entries in ``self`` equal to -1.
630
631
EXAMPLES::
632
633
sage: A = AlternatingSignMatrices(3)
634
sage: asm = A([[0,1,0],[1,0,0],[0,0,1]])
635
sage: asm.number_negative_ones()
636
0
637
sage: asm = A([[0,1,0],[1,-1,1],[0,1,0]])
638
sage: asm.number_negative_ones()
639
1
640
"""
641
a = self._matrix
642
return sum(1 for (i,j) in a.nonzero_positions() if a[i,j] == -1)
643
644
def is_permutation(self):
645
"""
646
Return ``True`` if ``self`` is a permutation matrix
647
and ``False`` otherwise.
648
649
EXAMPLES::
650
651
sage: A = AlternatingSignMatrices(3)
652
sage: asm = A([[0,1,0],[1,0,0],[0,0,1]])
653
sage: asm.is_permutation()
654
True
655
sage: asm = A([[0,1,0],[1,-1,1],[0,1,0]])
656
sage: asm.is_permutation()
657
False
658
"""
659
return self.number_negative_ones() == 0
660
661
def to_permutation(self):
662
"""
663
Return the corresponding permutation if ``self`` is a permutation
664
matrix.
665
666
EXAMPLES::
667
668
sage: A = AlternatingSignMatrices(3)
669
sage: asm = A([[0,1,0],[1,0,0],[0,0,1]])
670
sage: p = asm.to_permutation(); p
671
[2, 1, 3]
672
sage: parent(p)
673
Standard permutations
674
sage: asm = A([[0,1,0],[1,-1,1],[0,1,0]])
675
sage: asm.to_permutation()
676
Traceback (most recent call last):
677
...
678
ValueError: Not a permutation matrix
679
"""
680
if not self.is_permutation():
681
raise ValueError('Not a permutation matrix')
682
asm_matrix = self.to_matrix()
683
return Permutation([ j+1 for (i,j) in asm_matrix.nonzero_positions() ])
684
685
@combinatorial_map(name='to semistandard tableau')
686
def to_semistandard_tableau(self):
687
"""
688
Return the semistandard tableau corresponding the monotone triangle
689
corresponding to ``self``.
690
691
EXAMPLES::
692
693
sage: A = AlternatingSignMatrices(3)
694
sage: A([[0,0,1],[1,0,0],[0,1,0]]).to_semistandard_tableau()
695
[[1, 1, 3], [2, 3], [3]]
696
sage: t = A([[0,1,0],[1,-1,1],[0,1,0]]).to_semistandard_tableau(); t
697
[[1, 1, 2], [2, 3], [3]]
698
sage: parent(t)
699
Semistandard tableaux
700
"""
701
from sage.combinat.tableau import SemistandardTableau, SemistandardTableaux
702
mt = self.to_monotone_triangle()
703
ssyt = [[0]*(len(mt) - j) for j in range(len(mt))]
704
for i in range(len(mt)):
705
for j in range(len(mt[i])):
706
ssyt[i][j] = mt[j][-(i+1)]
707
return SemistandardTableau(ssyt)
708
709
def left_key(self):
710
r"""
711
Return the left key of the alternating sign matrix ``self``.
712
713
The left key of an alternating sign matrix was defined by Lascoux
714
in [LascouxPreprint]_ and is obtained by successively removing all the
715
`-1`'suntil what remains is a permutation matrix. This notion
716
corresponds to the notion of left key for semistandard tableaux. So
717
our algorithm proceeds as follows: we map ``self`` to its
718
corresponding monotone triangle, view that monotone triangle as a
719
semistandard tableaux, take its left key, and then map back through
720
monotone triangles to the permutation matrix which is the left key.
721
722
REFERENCES:
723
724
.. [Aval07] J.-C. Aval. *Keys and alternating sign matrices*.
725
Sem. Lothar. Combin. 59 (2007/10), Art. B59f, 13 pp.
726
727
.. [LascouxPreprint] A. Lascoux. *Chern and Yang through ice*.
728
Preprint.
729
730
EXAMPLES::
731
732
sage: A = AlternatingSignMatrices(3)
733
sage: A([[0,0,1],[1,0,0],[0,1,0]]).left_key()
734
[0 0 1]
735
[1 0 0]
736
[0 1 0]
737
sage: t = A([[0,1,0],[1,-1,1],[0,1,0]]).left_key(); t
738
[1 0 0]
739
[0 0 1]
740
[0 1 0]
741
sage: parent(t)
742
Alternating sign matrices of size 3
743
"""
744
from sage.combinat.tableau import SemistandardTableau, SemistandardTableaux
745
lkey = self.to_semistandard_tableau().left_key_tableau()
746
mt = [[0]*(len(lkey) - j) for j in range(len(lkey))]
747
for i in range(len(lkey)):
748
for j in range(len(lkey[i])):
749
mt[i][j] = lkey[len(lkey[i])-j-1][i]
750
A = AlternatingSignMatrices(len(lkey))
751
return A.from_monotone_triangle(mt)
752
753
@combinatorial_map(name='to left key permutation')
754
def left_key_as_permutation(self):
755
"""
756
Return the permutation of the left key of ``self``.
757
758
.. SEEALSO::
759
760
- :meth:`left_key()`
761
762
EXAMPLES::
763
764
sage: A = AlternatingSignMatrices(3)
765
sage: A([[0,0,1],[1,0,0],[0,1,0]]).left_key_as_permutation()
766
[3, 1, 2]
767
sage: t = A([[0,1,0],[1,-1,1],[0,1,0]]).left_key_as_permutation(); t
768
[1, 3, 2]
769
sage: parent(t)
770
Standard permutations
771
"""
772
return self.left_key().to_permutation()
773
774
class AlternatingSignMatrices(Parent, UniqueRepresentation):
775
r"""
776
Class of all `n \times n` alternating sign matrices.
777
778
An alternating sign matrix of size `n` is an `n \times n` matrix of `0`'s,
779
`1`'s and `-1`'s such that the sum of each row and column is `1` and the
780
non-zero entries in each row and column alternate in sign.
781
782
Alternating sign matrices of size `n` are in bijection with
783
:class:`monotone triangles <MonotoneTriangles>` with `n` rows.
784
785
INPUT:
786
787
- `n` -- an integer, the size of the matrices.
788
789
- ``use_monotone_triangle`` -- (Default: ``True``) If ``True``, the
790
generation of the matrices uses monotone triangles, else it will use the
791
earlier and now obsolete contre-tableaux implementation;
792
must be ``True`` to generate a lattice (with the ``lattice`` method)
793
794
EXAMPLES:
795
796
This will create an instance to manipulate the alternating sign
797
matrices of size 3::
798
799
sage: A = AlternatingSignMatrices(3)
800
sage: A
801
Alternating sign matrices of size 3
802
sage: A.cardinality()
803
7
804
805
Notably, this implementation allows to make a lattice of it::
806
807
sage: L = A.lattice()
808
sage: L
809
Finite lattice containing 7 elements
810
sage: L.category()
811
Category of facade finite lattice posets
812
"""
813
def __init__(self, n, use_monotone_triangles=True):
814
r"""
815
Initialize ``self``.
816
817
TESTS::
818
819
sage: A = AlternatingSignMatrices(4)
820
sage: TestSuite(A).run()
821
sage: A == AlternatingSignMatrices(4, use_monotone_triangles=False)
822
False
823
"""
824
self._n = n
825
self._matrix_space = MatrixSpace(ZZ, n)
826
self._umt = use_monotone_triangles
827
Parent.__init__(self, category=FiniteEnumeratedSets())
828
829
def _repr_(self):
830
r"""
831
Return a string representation of ``self``.
832
833
TESTS::
834
835
sage: A = AlternatingSignMatrices(4); A
836
Alternating sign matrices of size 4
837
"""
838
return "Alternating sign matrices of size %s" % self._n
839
840
def _repr_option(self, key):
841
"""
842
Metadata about the :meth:`_repr_` output.
843
844
See :meth:`sage.structure.parent._repr_option` for details.
845
846
EXAMPLES::
847
848
sage: A = AlternatingSignMatrices(3)
849
sage: A._repr_option('element_ascii_art')
850
True
851
"""
852
return self._matrix_space._repr_option(key)
853
854
def __contains__(self, asm):
855
"""
856
Check if ``asm`` is in ``self``.
857
858
TESTS::
859
860
sage: A = AlternatingSignMatrices(3)
861
sage: [[0,1,0],[1,0,0],[0,0,1]] in A
862
True
863
sage: [[0,1,0],[1,-1,1],[0,1,0]] in A
864
True
865
sage: [[0, 1],[1,0]] in A
866
False
867
sage: [[1,0,0,0],[0,1,0,0],[0,0,1,0],[0,0,0,1]] in A
868
False
869
sage: [[-1, 1, 1],[1,-1,1],[1,1,-1]] in A
870
False
871
"""
872
if isinstance(asm, AlternatingSignMatrix):
873
return asm._matrix.nrows() == self._n
874
try:
875
asm = self._matrix_space(asm)
876
except (TypeError, ValueError):
877
return False
878
for row in asm:
879
pos = False
880
for val in row:
881
if val > 0:
882
if pos:
883
return False
884
else:
885
pos = True
886
elif val < 0:
887
if pos:
888
pos = False
889
else:
890
return False
891
if not pos:
892
return False
893
if any(sum(row) != 1 for row in asm.columns()):
894
return False
895
return True
896
897
def _element_constructor_(self, asm):
898
"""
899
Construct an element of ``self``.
900
901
EXAMPLES::
902
903
sage: A = AlternatingSignMatrices(3)
904
sage: elt = A([[1, 0, 0],[0, 1, 0],[0, 0, 1]]); elt
905
[1 0 0]
906
[0 1 0]
907
[0 0 1]
908
sage: elt.parent() is A
909
True
910
sage: A([[3, 2, 1], [2, 1], [1]])
911
[1 0 0]
912
[0 1 0]
913
[0 0 1]
914
"""
915
if isinstance(asm, AlternatingSignMatrix):
916
if asm.parent() is self:
917
return asm
918
raise ValueError("Cannot convert between alternating sign matrices of different sizes")
919
if asm in MonotoneTriangles(self._n):
920
return self.from_monotone_triangle(asm)
921
return self.element_class(self, self._matrix_space(asm))
922
923
Element = AlternatingSignMatrix
924
925
def _an_element_(self):
926
"""
927
Return an element of ``self``.
928
"""
929
return self.element_class(self, self._matrix_space.identity_matrix())
930
931
def from_monotone_triangle(self, triangle):
932
r"""
933
Return an alternating sign matrix from a monotone triangle.
934
935
EXAMPLES::
936
937
sage: A = AlternatingSignMatrices(3)
938
sage: A.from_monotone_triangle([[3, 2, 1], [2, 1], [1]])
939
[1 0 0]
940
[0 1 0]
941
[0 0 1]
942
sage: A.from_monotone_triangle([[3, 2, 1], [3, 2], [3]])
943
[0 0 1]
944
[0 1 0]
945
[1 0 0]
946
"""
947
n = len(triangle)
948
if n != self._n:
949
raise ValueError("Incorrect size")
950
asm = []
951
952
prev = [0]*n
953
for line in reversed(triangle):
954
v = [1 if j+1 in reversed(line) else 0 for j in range(n)]
955
row = [a-b for (a, b) in zip(v, prev)]
956
asm.append(row)
957
prev = v
958
959
return self.element_class(self, self._matrix_space(asm))
960
961
def from_corner_sum(self, corner):
962
r"""
963
Return an alternating sign matrix from a corner sum matrix.
964
965
EXAMPLES::
966
967
sage: A = AlternatingSignMatrices(3)
968
sage: A.from_corner_sum(matrix([[0,0,0,0],[0,1,1,1],[0,1,2,2],[0,1,2,3]]))
969
[1 0 0]
970
[0 1 0]
971
[0 0 1]
972
sage: A.from_corner_sum(matrix([[0,0,0,0],[0,0,1,1],[0,1,1,2],[0,1,2,3]]))
973
[ 0 1 0]
974
[ 1 -1 1]
975
[ 0 1 0]
976
"""
977
asm_list=[]
978
n = len(list(corner)) - 1
979
for k in range(n):
980
asm_list.append([])
981
for i in range(n):
982
for j in range(n):
983
y = corner[i+1][j+1] \
984
- sum([sum([asm_list[i2][j2] for i2 in range(i)])
985
for j2 in range(j)]) \
986
- sum([asm_list[i2][j] for i2 in range(i)]) \
987
- sum([asm_list[i][j2] for j2 in range(j)])
988
asm_list[i].append(y)
989
return AlternatingSignMatrix(asm_list)
990
991
def from_height_function(self,height):
992
r"""
993
Return an alternating sign matrix from a height function.
994
995
EXAMPLES::
996
997
sage: A = AlternatingSignMatrices(3)
998
sage: A.from_height_function(matrix([[0,1,2,3],[1,2,1,2],[2,3,2,1],[3,2,1,0]]))
999
[0 0 1]
1000
[1 0 0]
1001
[0 1 0]
1002
sage: A.from_height_function(matrix([[0,1,2,3],[1,2,1,2],[2,1,2,1],[3,2,1,0]]))
1003
[ 0 1 0]
1004
[ 1 -1 1]
1005
[ 0 1 0]
1006
"""
1007
return self.from_corner_sum(matrix( [[((i+j-height[i][j])/int(2))
1008
for i in range(len(list(height)))]
1009
for j in range(len(list(height)))] ))
1010
1011
def size(self):
1012
r"""
1013
Return the size of the matrices in ``self``.
1014
1015
TESTS::
1016
1017
sage: A = AlternatingSignMatrices(4)
1018
sage: A.size()
1019
4
1020
"""
1021
return self._n
1022
1023
def cardinality(self):
1024
r"""
1025
Return the cardinality of ``self``.
1026
1027
The number of `n \times n` alternating sign matrices is equal to
1028
1029
.. MATH::
1030
1031
\prod_{k=0}^{n-1} \frac{(3k+1)!}{(n+k)!} = \frac{1! 4! 7! 10!
1032
\cdots (3n-2)!}{n! (n+1)! (n+2)! (n+3)! \cdots (2n-1)!}
1033
1034
EXAMPLES::
1035
1036
sage: [AlternatingSignMatrices(n).cardinality() for n in range(0, 11)]
1037
[1, 1, 2, 7, 42, 429, 7436, 218348, 10850216, 911835460, 129534272700]
1038
"""
1039
return Integer(prod( [ factorial(3*k+1)/factorial(self._n+k)
1040
for k in range(self._n)] ))
1041
1042
def matrix_space(self):
1043
"""
1044
Return the underlying matrix space.
1045
1046
EXAMPLES::
1047
1048
sage: A = AlternatingSignMatrices(3)
1049
sage: A.matrix_space()
1050
Full MatrixSpace of 3 by 3 dense matrices over Integer Ring
1051
"""
1052
return self._matrix_space
1053
1054
def __iter__(self):
1055
r"""
1056
Iterator on the alternating sign matrices of size `n`.
1057
1058
If defined using ``use_monotone_triangles``, this iterator
1059
will use the iteration on the monotone triangles. Else, it
1060
will use the iteration on contre-tableaux.
1061
1062
TESTS::
1063
1064
sage: A = AlternatingSignMatrices(4)
1065
sage: len(list(A))
1066
42
1067
"""
1068
if self._umt:
1069
for t in MonotoneTriangles(self._n):
1070
yield self.from_monotone_triangle(t)
1071
else:
1072
for c in ContreTableaux(self._n):
1073
yield from_contre_tableau(c)
1074
1075
def _lattice_initializer(self):
1076
r"""
1077
Return a 2-tuple to use in argument of ``LatticePoset``.
1078
1079
For more details about the cover relations, see
1080
``MonotoneTriangles``. Notice that the returned matrices are
1081
made immutable to ensure their hashability required by
1082
``LatticePoset``.
1083
1084
EXAMPLES:
1085
1086
Proof of the lattice property for alternating sign matrices of
1087
size 3::
1088
1089
sage: A = AlternatingSignMatrices(3)
1090
sage: P = Poset(A._lattice_initializer())
1091
sage: P.is_lattice()
1092
True
1093
"""
1094
assert(self._umt)
1095
(mts, rels) = MonotoneTriangles(self._n)._lattice_initializer()
1096
bij = dict((t, self.from_monotone_triangle(t)) for t in mts)
1097
asms, rels = bij.itervalues(), [(bij[a], bij[b]) for (a,b) in rels]
1098
return (asms, rels)
1099
1100
def cover_relations(self):
1101
r"""
1102
Iterate on the cover relations between the alternating sign
1103
matrices.
1104
1105
EXAMPLES::
1106
1107
sage: A = AlternatingSignMatrices(3)
1108
sage: for (a,b) in A.cover_relations():
1109
....: eval('a, b')
1110
(
1111
[1 0 0] [0 1 0]
1112
[0 1 0] [1 0 0]
1113
[0 0 1], [0 0 1]
1114
)
1115
(
1116
[1 0 0] [1 0 0]
1117
[0 1 0] [0 0 1]
1118
[0 0 1], [0 1 0]
1119
)
1120
(
1121
[0 1 0] [ 0 1 0]
1122
[1 0 0] [ 1 -1 1]
1123
[0 0 1], [ 0 1 0]
1124
)
1125
(
1126
[1 0 0] [ 0 1 0]
1127
[0 0 1] [ 1 -1 1]
1128
[0 1 0], [ 0 1 0]
1129
)
1130
(
1131
[ 0 1 0] [0 0 1]
1132
[ 1 -1 1] [1 0 0]
1133
[ 0 1 0], [0 1 0]
1134
)
1135
(
1136
[ 0 1 0] [0 1 0]
1137
[ 1 -1 1] [0 0 1]
1138
[ 0 1 0], [1 0 0]
1139
)
1140
(
1141
[0 0 1] [0 0 1]
1142
[1 0 0] [0 1 0]
1143
[0 1 0], [1 0 0]
1144
)
1145
(
1146
[0 1 0] [0 0 1]
1147
[0 0 1] [0 1 0]
1148
[1 0 0], [1 0 0]
1149
)
1150
1151
"""
1152
(_, rels) = self._lattice_initializer()
1153
return (_ for _ in rels)
1154
1155
def lattice(self):
1156
r"""
1157
Return the lattice of the alternating sign matrices of size
1158
`n`, created by ``LatticePoset``.
1159
1160
EXAMPLES::
1161
1162
sage: A = AlternatingSignMatrices(3)
1163
sage: L = A.lattice()
1164
sage: L
1165
Finite lattice containing 7 elements
1166
1167
"""
1168
return LatticePoset(self._lattice_initializer(), cover_relations=True)
1169
1170
1171
class MonotoneTriangles(GelfandTsetlinPatternsTopRow):
1172
r"""
1173
Monotone triangles with `n` rows.
1174
1175
A monotone triangle is a number triangle `(a_{i,j})_{1 \leq i \leq
1176
n , 1 \leq j \leq i}` on `\{1, \dots, n\}` such that:
1177
1178
- `a_{i,j} < a_{i,j+1}`
1179
1180
- `a_{i+1,j} < a_{i,j} \leq a_{i+1,j+1}`
1181
1182
This notably requires that the bottom column is ``[1,...,n]``.
1183
1184
Alternatively a monotone triangle is a strict Gelfand-Tsetlin pattern with
1185
top row `(n, \ldots, 2, 1)`.
1186
1187
INPUT:
1188
1189
- ``n`` -- The number of rows in the monotone triangles
1190
1191
EXAMPLES:
1192
1193
This represents the monotone triangles with base ``[3,2,1]``::
1194
1195
sage: M = MonotoneTriangles(3)
1196
sage: M
1197
Monotone triangles with 3 rows
1198
sage: M.cardinality()
1199
7
1200
1201
The monotone triangles are a lattice::
1202
1203
sage: M.lattice()
1204
Finite lattice containing 7 elements
1205
1206
Monotone triangles can be converted to alternating sign matrices
1207
and back::
1208
1209
sage: M = MonotoneTriangles(5)
1210
sage: A = AlternatingSignMatrices(5)
1211
sage: all(A.from_monotone_triangle(m).to_monotone_triangle() == m for m in M)
1212
True
1213
"""
1214
def __init__(self, n):
1215
r"""
1216
Initialize ``self``.
1217
1218
TESTS::
1219
1220
sage: M = MonotoneTriangles(4)
1221
sage: TestSuite(M).run()
1222
sage: M2 = MonotoneTriangles(int(4))
1223
sage: M is M2
1224
True
1225
"""
1226
GelfandTsetlinPatternsTopRow.__init__(self, tuple(reversed(range(1, n+1))), True)
1227
1228
def _repr_(self):
1229
r"""
1230
String representation.
1231
1232
TESTS::
1233
1234
sage: M = MonotoneTriangles(4)
1235
sage: M
1236
Monotone triangles with 4 rows
1237
"""
1238
return "Monotone triangles with %s rows" % self._n
1239
1240
def cardinality(self):
1241
r"""
1242
Cardinality of ``self``.
1243
1244
The number of monotone triangles with `n` rows is equal to
1245
1246
.. MATH::
1247
1248
\prod_{k=0}^{n-1} \frac{(3k+1)!}{(n+k)!} = \frac{1! 4! 7! 10!
1249
\cdots (3n-2)!}{n! (n+1)! (n+2)! (n+3)! \cdots (2n-1)!}
1250
1251
EXAMPLES::
1252
1253
sage: M = MonotoneTriangles(4)
1254
sage: M.cardinality()
1255
42
1256
"""
1257
return Integer(prod( [ factorial(3*k+1)/factorial(self._n+k)
1258
for k in range(self._n)] ))
1259
1260
def _lattice_initializer(self):
1261
r"""
1262
Return a 2-tuple to use in argument of ``LatticePoset``.
1263
1264
This couple is composed by the set of the monotone triangles
1265
with `n` rows and the cover relations. Specializing this
1266
function allows to generate the monotone triangles just once,
1267
and so to speed up the computation in comparison of
1268
``(list(self), self.cover_relations())``. Notice that the
1269
function also switch the representation of monotone triangles
1270
from list of list to tuple of tuple in order to make them
1271
hashable (required to make a poset with them).
1272
1273
EXAMPLES::
1274
1275
sage: M = MonotoneTriangles(3)
1276
sage: P = Poset(M._lattice_initializer())
1277
sage: P.is_lattice()
1278
True
1279
"""
1280
# get a list of the elements and switch to a tuple
1281
# representation
1282
set_ = list(self)
1283
set_ = map(lambda x: tuple(map(tuple, x)), set_)
1284
return (set_, [(a,b) for a in set_ for b in set_ if _is_a_cover(a,b)])
1285
1286
def cover_relations(self):
1287
r"""
1288
Iterate on the cover relations in the set of monotone triangles
1289
with `n` rows.
1290
1291
EXAMPLES::
1292
1293
sage: M = MonotoneTriangles(3)
1294
sage: for (a,b) in M.cover_relations():
1295
....: eval('a, b')
1296
([[3, 2, 1], [2, 1], [1]], [[3, 2, 1], [2, 1], [2]])
1297
([[3, 2, 1], [2, 1], [1]], [[3, 2, 1], [3, 1], [1]])
1298
([[3, 2, 1], [2, 1], [2]], [[3, 2, 1], [3, 1], [2]])
1299
([[3, 2, 1], [3, 1], [1]], [[3, 2, 1], [3, 1], [2]])
1300
([[3, 2, 1], [3, 1], [2]], [[3, 2, 1], [3, 1], [3]])
1301
([[3, 2, 1], [3, 1], [2]], [[3, 2, 1], [3, 2], [2]])
1302
([[3, 2, 1], [3, 1], [3]], [[3, 2, 1], [3, 2], [3]])
1303
([[3, 2, 1], [3, 2], [2]], [[3, 2, 1], [3, 2], [3]])
1304
"""
1305
set_ = list(self)
1306
return ((a,b) for a in set_ for b in set_ if _is_a_cover(a,b))
1307
1308
def lattice(self):
1309
r"""
1310
Return the lattice of the monotone triangles with `n` rows.
1311
1312
EXAMPLES::
1313
1314
sage: M = MonotoneTriangles(3)
1315
sage: P = M.lattice()
1316
sage: P
1317
Finite lattice containing 7 elements
1318
1319
"""
1320
return LatticePoset(self._lattice_initializer(), cover_relations=True)
1321
1322
def _is_a_cover(mt0, mt1):
1323
r"""
1324
Define the cover relations.
1325
1326
Return ``True`` if and only if the second argument is a cover of
1327
the first one.
1328
1329
EXAMPLES::
1330
1331
sage: import sage.combinat.alternating_sign_matrix as asm
1332
sage: asm._is_a_cover([[1,2,3],[1,2],[1]], [[1,2,3],[1,3],[1]])
1333
True
1334
sage: asm._is_a_cover([[1,2,3],[1,3],[2]], [[1,2,3],[1,2],[1]])
1335
False
1336
"""
1337
diffs = 0
1338
for (a,b) in itertools.izip(flatten(mt0), flatten(mt1)):
1339
if a != b:
1340
if a+1 == b:
1341
diffs += 1
1342
else:
1343
return False
1344
if diffs > 1:
1345
return False
1346
return diffs == 1
1347
1348
# Deprecated methods
1349
1350
def to_monotone_triangle(matrix):
1351
"""
1352
Deprecated method, use :meth:`AlternatingSignMatrix.to_monotone_triangle()`
1353
instead.
1354
1355
EXAMPLES::
1356
1357
sage: sage.combinat.alternating_sign_matrix.to_monotone_triangle([[0,1],[1,0]])
1358
doctest:...: DeprecationWarning: to_monotone_triangle() is deprecated. Use AlternatingSignMatrix.to_monotone_triangle() instead
1359
See http://trac.sagemath.org/14301 for details.
1360
[[2, 1], [2]]
1361
"""
1362
from sage.misc.superseded import deprecation
1363
deprecation(14301,'to_monotone_triangle() is deprecated. Use AlternatingSignMatrix.to_monotone_triangle() instead')
1364
return AlternatingSignMatrix(matrix).to_monotone_triangle()
1365
1366
def from_monotone_triangle(triangle):
1367
"""
1368
Deprecated method, use
1369
:meth:`AlternatingSignMatrices.from_monotone_triangle()` instead.
1370
1371
EXAMPLES::
1372
1373
sage: sage.combinat.alternating_sign_matrix.from_monotone_triangle([[1, 2], [2]])
1374
doctest:...: DeprecationWarning: from_monotone_triangle() is deprecated. Use AlternatingSignMatrix.from_monotone_triangle() instead
1375
See http://trac.sagemath.org/14301 for details.
1376
[0 1]
1377
[1 0]
1378
"""
1379
from sage.misc.superseded import deprecation
1380
deprecation(14301,'from_monotone_triangle() is deprecated. Use AlternatingSignMatrix.from_monotone_triangle() instead')
1381
return AlternatingSignMatrices(len(triangle)).from_monotone_triangle(triangle)
1382
1383
# For old pickles
1384
def AlternatingSignMatrices_n(n):
1385
"""
1386
For old pickles of ``AlternatingSignMatrices_n``.
1387
1388
EXAMPLES::
1389
1390
sage: sage.combinat.alternating_sign_matrix.AlternatingSignMatrices_n(3)
1391
doctest:...: DeprecationWarning: this class is deprecated. Use sage.combinat.alternating_sign_matrix.AlternatingSignMatrices instead
1392
See http://trac.sagemath.org/14301 for details.
1393
Alternating sign matrices of size 3
1394
"""
1395
from sage.misc.superseded import deprecation
1396
deprecation(14301,'this class is deprecated. Use sage.combinat.alternating_sign_matrix.AlternatingSignMatrices instead')
1397
return AlternatingSignMatrices(n)
1398
1399
def MonotoneTriangles_n(n):
1400
"""
1401
For old pickles of ``MonotoneTriangles_n``.
1402
1403
EXAMPLES::
1404
1405
sage: sage.combinat.alternating_sign_matrix.MonotoneTriangles_n(3)
1406
doctest:...: DeprecationWarning: this class is deprecated. Use sage.combinat.alternating_sign_matrix.MonotoneTriangles instead
1407
See http://trac.sagemath.org/14301 for details.
1408
Monotone triangles with 3 rows
1409
"""
1410
from sage.misc.superseded import deprecation
1411
deprecation(14301,'this class is deprecated. Use sage.combinat.alternating_sign_matrix.MonotoneTriangles instead')
1412
return MonotoneTriangles(n)
1413
1414
from sage.structure.sage_object import register_unpickle_override
1415
register_unpickle_override('sage.combinat.alternating_sign_matrix', 'AlternatingSignMatrices_n', AlternatingSignMatrices)
1416
register_unpickle_override('sage.combinat.alternating_sign_matrix', 'MonotoneTriangles_n', MonotoneTriangles)
1417
register_unpickle_override('sage.combinat.alternating_sign_matrix', 'MonotoneTriangles_n', MonotoneTriangles_n)
1418
1419
# Here are the previous implementations of the combinatorial structure
1420
# of the alternating sign matrices. Please, consider it obsolete and
1421
# tend to use the monotone triangles instead.
1422
1423
def from_contre_tableau(comps):
1424
r"""
1425
Returns an alternating sign matrix from a contre-tableau.
1426
1427
EXAMPLES::
1428
1429
sage: import sage.combinat.alternating_sign_matrix as asm
1430
sage: asm.from_contre_tableau([[1, 2, 3], [1, 2], [1]])
1431
doctest:...: DeprecationWarning: You can use from_monotone_triangle instead.
1432
See http://trac.sagemath.org/12930 for details.
1433
[0 0 1]
1434
[0 1 0]
1435
[1 0 0]
1436
sage: asm.from_contre_tableau([[1, 2, 3], [2, 3], [3]])
1437
[1 0 0]
1438
[0 1 0]
1439
[0 0 1]
1440
"""
1441
from sage.misc.superseded import deprecation
1442
deprecation(12930, 'You can use from_monotone_triangle instead.')
1443
n = len(comps)
1444
MS = MatrixSpace(ZZ, n)
1445
M = [ [0 for _ in range(n)] for _ in range(n) ]
1446
1447
previous_set = Set([])
1448
1449
for col in range(n-1, -1, -1):
1450
s = Set( comps[col] )
1451
for x in s - previous_set:
1452
M[x-1][col] = 1
1453
1454
for x in previous_set - s:
1455
M[x-1][col] = -1
1456
1457
previous_set = s
1458
1459
return MS(M)
1460
1461
1462
class ContreTableaux(Parent):
1463
"""
1464
Factory class for the combinatorial class of contre tableaux of size `n`.
1465
1466
EXAMPLES::
1467
1468
sage: ct4 = ContreTableaux(4); ct4
1469
Contre tableaux of size 4
1470
sage: ct4.cardinality()
1471
42
1472
"""
1473
__metaclass__ = ClasscallMetaclass
1474
@staticmethod
1475
def __classcall_private__(cls, n, **kwds):
1476
r"""
1477
Factory pattern.
1478
1479
Check properties on arguments, then call the appropriate class.
1480
1481
EXAMPLES::
1482
1483
sage: C = ContreTableaux(4)
1484
sage: type(C)
1485
<class 'sage.combinat.alternating_sign_matrix.ContreTableaux_n'>
1486
1487
"""
1488
assert(isinstance(n, (int, Integer)))
1489
return ContreTableaux_n(n, **kwds)
1490
1491
1492
class ContreTableaux_n(ContreTableaux):
1493
def __init__(self, n):
1494
"""
1495
TESTS::
1496
1497
sage: ct2 = ContreTableaux(2); ct2
1498
Contre tableaux of size 2
1499
sage: ct2 == loads(dumps(ct2))
1500
True
1501
"""
1502
self.n = n
1503
1504
def __repr__(self):
1505
"""
1506
TESTS::
1507
1508
sage: repr(ContreTableaux(2))
1509
'Contre tableaux of size 2'
1510
"""
1511
return "Contre tableaux of size %s"%self.n
1512
1513
def __eq__(self, other):
1514
"""
1515
TESTS::
1516
1517
sage: C = ContreTableaux(4)
1518
sage: C == loads(dumps(C))
1519
True
1520
1521
"""
1522
return self.n == other.n
1523
1524
def cardinality(self):
1525
"""
1526
EXAMPLES::
1527
1528
sage: [ ContreTableaux(n).cardinality() for n in range(0, 11)]
1529
[1, 1, 2, 7, 42, 429, 7436, 218348, 10850216, 911835460, 129534272700]
1530
"""
1531
return prod( [ factorial(3*k+1)/factorial(self.n+k) for k in range(self.n)] )
1532
1533
def _iterator_rec(self, i):
1534
"""
1535
EXAMPLES::
1536
1537
sage: c = ContreTableaux(2)
1538
sage: list(c._iterator_rec(0))
1539
[[]]
1540
sage: list(c._iterator_rec(1))
1541
[[[1, 2]]]
1542
sage: list(c._iterator_rec(2))
1543
[[[1, 2], [1]], [[1, 2], [2]]]
1544
"""
1545
if i == 0:
1546
yield []
1547
elif i == 1:
1548
yield [range(1, self.n+1)]
1549
else:
1550
for columns in self._iterator_rec(i-1):
1551
previous_column = columns[-1]
1552
for column in _next_column_iterator(previous_column, len(previous_column)-1):
1553
yield columns + [ column ]
1554
1555
def __iter__(self):
1556
"""
1557
EXAMPLES::
1558
1559
sage: list(ContreTableaux(0))
1560
[[]]
1561
sage: list(ContreTableaux(1))
1562
[[[1]]]
1563
sage: list(ContreTableaux(2))
1564
[[[1, 2], [1]], [[1, 2], [2]]]
1565
sage: list(ContreTableaux(3))
1566
[[[1, 2, 3], [1, 2], [1]],
1567
[[1, 2, 3], [1, 2], [2]],
1568
[[1, 2, 3], [1, 3], [1]],
1569
[[1, 2, 3], [1, 3], [2]],
1570
[[1, 2, 3], [1, 3], [3]],
1571
[[1, 2, 3], [2, 3], [2]],
1572
[[1, 2, 3], [2, 3], [3]]]
1573
"""
1574
1575
for z in self._iterator_rec(self.n):
1576
yield z
1577
1578
1579
def _next_column_iterator(previous_column, height, i = None):
1580
"""
1581
Returns a generator for all columns of height height properly
1582
filled from row 1 to ``i``
1583
1584
EXAMPLES::
1585
1586
sage: import sage.combinat.alternating_sign_matrix as asm
1587
sage: list(asm._next_column_iterator([1], 0))
1588
[[]]
1589
sage: list(asm._next_column_iterator([1,5],1))
1590
[[1], [2], [3], [4], [5]]
1591
sage: list(asm._next_column_iterator([1,4,5],2))
1592
[[1, 4], [1, 5], [2, 4], [2, 5], [3, 4], [3, 5], [4, 5]]
1593
"""
1594
if i is None:
1595
i = height
1596
if i == 0:
1597
yield [-1]*height
1598
else:
1599
for column in _next_column_iterator(previous_column, height, i-1):
1600
min_value = previous_column[i-1]
1601
if i > 1:
1602
min_value = max(min_value, column[i-2]+1)
1603
for value in range(min_value, previous_column[i]+1):
1604
c = copy.copy(column)
1605
c[i-1] = value
1606
yield c
1607
1608
1609
def _previous_column_iterator(column, height, max_value):
1610
"""
1611
EXAMPLES::
1612
1613
sage: import sage.combinat.alternating_sign_matrix as asm
1614
sage: list(asm._previous_column_iterator([2,3], 3, 4))
1615
[[1, 2, 3], [1, 2, 4], [1, 3, 4], [2, 3, 4]]
1616
"""
1617
new_column = [1] + column + [ max_value ] * (height - len(column))
1618
return _next_column_iterator(new_column, height)
1619
1620
1621
class TruncatedStaircases(Parent):
1622
"""
1623
Factory class for the combinatorial class of truncated staircases
1624
of size ``n`` with last column ``last_column``.
1625
1626
EXAMPLES::
1627
1628
sage: t4 = TruncatedStaircases(4, [2,3]); t4
1629
Truncated staircases of size 4 with last column [2, 3]
1630
sage: t4.cardinality()
1631
4
1632
"""
1633
__metaclass__ = ClasscallMetaclass
1634
@staticmethod
1635
def __classcall_private__(cls, n, last_column, **kwds):
1636
r"""
1637
Factory pattern.
1638
1639
Check properties on arguments, then call the appropriate class.
1640
1641
TESTS::
1642
1643
sage: T = TruncatedStaircases(4, [2,3])
1644
sage: type(T)
1645
<class 'sage.combinat.alternating_sign_matrix.TruncatedStaircases_nlastcolumn'>
1646
1647
"""
1648
assert(isinstance(n, (int, Integer)))
1649
return TruncatedStaircases_nlastcolumn(n, last_column, **kwds)
1650
1651
1652
class TruncatedStaircases_nlastcolumn(TruncatedStaircases):
1653
def __init__(self, n, last_column):
1654
"""
1655
TESTS::
1656
1657
sage: t4 = TruncatedStaircases(4, [2,3]); t4
1658
Truncated staircases of size 4 with last column [2, 3]
1659
sage: t4 == loads(dumps(t4))
1660
True
1661
"""
1662
self.n = n
1663
self.last_column = last_column
1664
1665
def __repr__(self):
1666
"""
1667
TESTS::
1668
1669
sage: repr(TruncatedStaircases(4, [2,3]))
1670
'Truncated staircases of size 4 with last column [2, 3]'
1671
"""
1672
return "Truncated staircases of size %s with last column %s"%(self.n, self.last_column)
1673
1674
def _iterator_rec(self, i):
1675
"""
1676
EXAMPLES::
1677
1678
sage: t = TruncatedStaircases(3, [2,3])
1679
sage: list(t._iterator_rec(1))
1680
[]
1681
sage: list(t._iterator_rec(2))
1682
[[[2, 3]]]
1683
sage: list(t._iterator_rec(3))
1684
[[[1, 2, 3], [2, 3]]]
1685
"""
1686
if i < len(self.last_column):
1687
return
1688
elif i == len(self.last_column):
1689
yield [self.last_column]
1690
else:
1691
for columns in self._iterator_rec(i-1):
1692
previous_column = columns[0]
1693
for column in _previous_column_iterator(previous_column, len(previous_column)+1, self.n):
1694
yield [column] + columns
1695
1696
def __iter__(self):
1697
"""
1698
EXAMPLES:::
1699
1700
sage: list(TruncatedStaircases(4, [2,3]))
1701
[[[4, 3, 2, 1], [3, 2, 1], [3, 2]], [[4, 3, 2, 1], [4, 2, 1], [3, 2]], [[4, 3, 2, 1], [4, 3, 1], [3, 2]], [[4, 3, 2, 1], [4, 3, 2], [3, 2]]]
1702
"""
1703
for z in self._iterator_rec(self.n):
1704
yield map(lambda x: list(reversed(x)), z)
1705
1706
def __eq__(self, other):
1707
r"""
1708
TESTS::
1709
1710
sage: T = TruncatedStaircases(4, [2,3])
1711
sage: T == loads(dumps(T))
1712
True
1713
"""
1714
return ((self.n == other.n) and
1715
(self.last_column == other.last_column))
1716
1717
def cardinality(self):
1718
r"""
1719
EXAMPLES::
1720
1721
sage: T = TruncatedStaircases(4, [2,3])
1722
sage: T.cardinality()
1723
4
1724
"""
1725
c = 0
1726
for _ in self:
1727
c += 1
1728
return c
1729
1730
def nw_corner_sum(M,i,j):
1731
r"""
1732
Return the sum of entries to the northwest of `(i,j)` in matrix.
1733
1734
EXAMPLES::
1735
1736
sage: from sage.combinat.alternating_sign_matrix import nw_corner_sum
1737
sage: A = matrix.ones(3,3)
1738
sage: nw_corner_sum(A,2,2)
1739
4
1740
"""
1741
if i >= 0 and j >= 0:
1742
return sum([sum([M[i2][j2] for j2 in range(j)]) for i2 in range(i)])
1743
return 0
1744
1745
1746