Book a Demo!
CoCalc Logo Icon
StoreFeaturesDocsShareSupportNewsAboutPoliciesSign UpSign In
sagemath
GitHub Repository: sagemath/sagesmc
Path: blob/master/src/sage/combinat/cluster_algebra_quiver/cluster_seed.py
8818 views
1
r"""
2
ClusterSeed
3
4
A *cluster seed* is a pair `(B,\mathbf{x})` with `B` being a *skew-symmetrizable* `(n+m \times n)` *-matrix*
5
and with `\mathbf{x}` being an `n`-tuple of *independent elements* in the field of rational functions in `n` variables.
6
7
For the compendium on the cluster algebra and quiver package see
8
:arxiv:`1102.4844`.
9
10
AUTHORS:
11
12
- Gregg Musiker
13
- Christian Stump
14
15
.. seealso:: For mutation types of cluster seeds, see :meth:`sage.combinat.cluster_algebra_quiver.quiver_mutation_type.QuiverMutationType`. Cluster seeds are closely related to :meth:`sage.combinat.cluster_algebra_quiver.quiver.ClusterQuiver`.
16
"""
17
18
#*****************************************************************************
19
# Copyright (C) 2011 Gregg Musiker <[email protected]>
20
# Christian Stump <[email protected]>
21
#
22
# Distributed under the terms of the GNU General Public License (GPL)
23
# http://www.gnu.org/licenses/
24
#*****************************************************************************
25
import time
26
from sage.structure.sage_object import SageObject
27
from copy import copy
28
from sage.rings.all import QQ, infinity
29
from sage.rings.all import FractionField, PolynomialRing
30
from sage.rings.fraction_field_element import FractionFieldElement
31
from sage.sets.all import Set
32
from sage.combinat.cluster_algebra_quiver.quiver_mutation_type import QuiverMutationType_Irreducible, QuiverMutationType_Reducible
33
from sage.combinat.cluster_algebra_quiver.mutation_type import is_mutation_finite
34
35
class ClusterSeed(SageObject):
36
r"""
37
The *cluster seed* associated to an *exchange matrix*.
38
39
INPUT:
40
41
- ``data`` -- can be any of the following::
42
43
* QuiverMutationType
44
* str - a string representing a QuiverMutationType or a common quiver type (see Examples)
45
* ClusterQuiver
46
* Matrix - a skew-symmetrizable matrix
47
* DiGraph - must be the input data for a quiver
48
* List of edges - must be the edge list of a digraph for a quiver
49
50
EXAMPLES::
51
52
sage: S = ClusterSeed(['A',5]); S
53
A seed for a cluster algebra of rank 5 of type ['A', 5]
54
55
sage: S = ClusterSeed(['A',[2,5],1]); S
56
A seed for a cluster algebra of rank 7 of type ['A', [2, 5], 1]
57
58
sage: T = ClusterSeed( S ); T
59
A seed for a cluster algebra of rank 7 of type ['A', [2, 5], 1]
60
61
sage: T = ClusterSeed( S._M ); T
62
A seed for a cluster algebra of rank 7
63
64
sage: T = ClusterSeed( S.quiver()._digraph ); T
65
A seed for a cluster algebra of rank 7
66
67
sage: T = ClusterSeed( S.quiver()._digraph.edges() ); T
68
A seed for a cluster algebra of rank 7
69
70
sage: S = ClusterSeed(['B',2]); S
71
A seed for a cluster algebra of rank 2 of type ['B', 2]
72
73
sage: S = ClusterSeed(['C',2]); S
74
A seed for a cluster algebra of rank 2 of type ['B', 2]
75
76
sage: S = ClusterSeed(['A', [5,0],1]); S
77
A seed for a cluster algebra of rank 5 of type ['D', 5]
78
79
sage: S = ClusterSeed(['GR',[3,7]]); S
80
A seed for a cluster algebra of rank 6 of type ['E', 6]
81
82
sage: S = ClusterSeed(['F', 4, [2,1]]); S
83
A seed for a cluster algebra of rank 6 of type ['F', 4, [1, 2]]
84
"""
85
def __init__(self, data, frozen=None, is_principal=None):
86
r"""
87
TESTS::
88
89
sage: S = ClusterSeed(['A',4])
90
sage: TestSuite(S).run()
91
"""
92
from quiver import ClusterQuiver
93
94
# constructs a cluster seed from a cluster seed
95
if type(data) is ClusterSeed:
96
if frozen:
97
print "The input \'frozen\' is ignored"
98
self._M = copy( data._M )
99
self._cluster = copy(data._cluster)
100
self._n = data._n
101
self._m = data._m
102
self._R = data._R
103
self._quiver = ClusterQuiver( data._quiver ) if data._quiver else None
104
self._mutation_type = copy( data._mutation_type )
105
self._description = copy( data._description )
106
self._is_principal = data._is_principal
107
108
# constructs a cluster seed from a quiver
109
elif type(data) is ClusterQuiver:
110
if frozen:
111
print "The input \'frozen\' is ignored"
112
113
quiver = ClusterQuiver( data )
114
self._M = copy(quiver._M)
115
self._n = quiver._n
116
self._m = quiver._m
117
self._quiver = quiver
118
self._mutation_type = quiver._mutation_type
119
self._description = 'A seed for a cluster algebra of rank %d' %(self._n)
120
self._R = FractionField(PolynomialRing(QQ,['x%s'%i for i in range(0,self._n)]+['y%s'%i for i in range(0,self._m)]))
121
self._cluster = list(self._R.gens())
122
self._is_principal = None
123
124
# in all other cases, we construct the corresponding ClusterQuiver first
125
else:
126
quiver = ClusterQuiver( data, frozen=frozen )
127
self.__init__( quiver )
128
129
if is_principal != None:
130
self._is_principal = is_principal
131
132
def __eq__(self, other):
133
r"""
134
Returns True iff ``self`` represent the same cluster seed as ``other``.
135
136
EXAMPLES::
137
138
sage: S = ClusterSeed(['A',5])
139
sage: T = S.mutate( 2, inplace=False )
140
sage: S.__eq__( T )
141
False
142
143
sage: T.mutate( 2 )
144
sage: S.__eq__( T )
145
True
146
"""
147
return type( other ) is ClusterSeed and self._M == other._M and self._cluster == other._cluster
148
149
def _repr_(self):
150
r"""
151
Returns the description of ``self``.
152
153
EXAMPLES::
154
155
sage: S = ClusterSeed(['A',5])
156
sage: S._repr_()
157
"A seed for a cluster algebra of rank 5 of type ['A', 5]"
158
159
sage: S=ClusterSeed(['B',2])
160
sage: T=S.principal_extension()
161
sage: T._repr_()
162
"A seed for a cluster algebra of rank 2 of type ['B', 2] with principal coefficients"
163
"""
164
name = self._description
165
if self._mutation_type:
166
if type( self._mutation_type ) in [QuiverMutationType_Irreducible,QuiverMutationType_Reducible]:
167
name += ' of type ' + str(self._mutation_type)
168
# the following case allows description of 'undetermined finite mutation type'
169
else:
170
name += ' of ' + self._mutation_type
171
if self._is_principal:
172
name += ' with principal coefficients'
173
elif self._m == 1:
174
name += ' with %s frozen variable'%self._m
175
elif self._m > 1:
176
name += ' with %s frozen variables'%self._m
177
return name
178
179
def plot(self, circular=False, mark=None, save_pos=False):
180
r"""
181
Returns the plot of the quiver of ``self``.
182
183
INPUT:
184
185
- ``circular`` -- (default:False) if True, the circular plot is chosen, otherwise >>spring<< is used.
186
- ``mark`` -- (default: None) if set to i, the vertex i is highlighted.
187
- ``save_pos`` -- (default:False) if True, the positions of the vertices are saved.
188
189
EXAMPLES::
190
191
sage: S = ClusterSeed(['A',5])
192
sage: pl = S.plot()
193
sage: pl = S.plot(circular=True)
194
"""
195
return self.quiver().plot(circular=circular,mark=mark,save_pos=save_pos)
196
197
def show(self, fig_size=1, circular=False, mark=None, save_pos=False):
198
r"""
199
Shows the plot of the quiver of ``self``.
200
201
INPUT:
202
203
- ``fig_size`` -- (default: 1) factor by which the size of the plot is multiplied.
204
- ``circular`` -- (default: False) if True, the circular plot is chosen, otherwise >>spring<< is used.
205
- ``mark`` -- (default: None) if set to i, the vertex i is highlighted.
206
- ``save_pos`` -- (default:False) if True, the positions of the vertices are saved.
207
208
TESTS::
209
210
sage: S = ClusterSeed(['A',5])
211
sage: S.show() # long time
212
"""
213
self.quiver().show(fig_size=fig_size, circular=circular,mark=mark,save_pos=save_pos)
214
215
def interact(self, fig_size=1, circular=True):
216
r"""
217
Only in *notebook mode*. Starts an interactive window for cluster seed mutations.
218
219
INPUT:
220
221
- ``fig_size`` -- (default: 1) factor by which the size of the plot is multiplied.
222
- ``circular`` -- (default: True) if True, the circular plot is chosen, otherwise >>spring<< is used.
223
224
TESTS::
225
226
sage: S = ClusterSeed(['A',4])
227
sage: S.interact() # long time
228
'The interactive mode only runs in the Sage notebook.'
229
"""
230
from sage.plot.plot import EMBEDDED_MODE
231
from sagenb.notebook.interact import interact, selector
232
from sage.misc.all import html,latex
233
from sage.all import var
234
235
if not EMBEDDED_MODE:
236
return "The interactive mode only runs in the Sage notebook."
237
else:
238
seq = []
239
sft = [True]
240
sss = [True]
241
ssv = [True]
242
ssm = [True]
243
ssl = [True]
244
@interact
245
def player(k=selector(values=range(self._n),nrows = 1,label='Mutate at: '), show_seq=("Mutation sequence:", True), show_vars=("Cluster variables:", True), show_matrix=("B-Matrix:", True), show_lastmutation=("Show last mutation:", True) ):
246
ft,ss,sv,sm,sl = sft.pop(), sss.pop(), ssv.pop(), ssm.pop(), ssl.pop()
247
if ft:
248
self.show(fig_size=fig_size, circular=circular)
249
elif show_seq is not ss or show_vars is not sv or show_matrix is not sm or show_lastmutation is not sl:
250
if seq and show_lastmutation:
251
self.show(fig_size=fig_size, circular=circular, mark=seq[len(seq)-1])
252
else:
253
self.show(fig_size=fig_size, circular=circular )
254
else:
255
self.mutate(k)
256
seq.append(k)
257
if not show_lastmutation:
258
self.show(fig_size=fig_size, circular=circular)
259
else:
260
self.show(fig_size=fig_size, circular=circular,mark=k)
261
sft.append(False)
262
sss.append(show_seq)
263
ssv.append(show_vars)
264
ssm.append(show_matrix)
265
ssl.append(show_lastmutation)
266
if show_seq: html( "Mutation sequence: $" + str( [ seq[i] for i in xrange(len(seq)) ] ).strip('[]') + "$" )
267
if show_vars:
268
html( "Cluster variables:" )
269
table = "$\\begin{align*}\n"
270
for i in xrange(self._n):
271
v = var('v%s'%i)
272
table += "\t" + latex( v ) + " &= " + latex( self._cluster[i] ) + "\\\\ \\\\\n"
273
table += "\\end{align*}$"
274
html( "$ $" )
275
html( table )
276
html( "$ $" )
277
if show_matrix:
278
html( "B-Matrix:" )
279
m = self._M
280
#m = matrix(range(1,self._n+1),sparse=True).stack(m)
281
m = latex(m)
282
m = m.split('(')[1].split('\\right')[0]
283
html( "$ $" )
284
html( "$\\begin{align*} " + m + "\\end{align*}$" )
285
#html( "$" + m + "$" )
286
html( "$ $" )
287
288
def save_image(self, filename, circular=False, mark=None, save_pos=False):
289
r"""
290
Saves the plot of the underlying digraph of the quiver of ``self``.
291
292
INPUT:
293
294
- ``filename`` -- the filename the image is saved to.
295
- ``circular`` -- (default: False) if True, the circular plot is chosen, otherwise >>spring<< is used.
296
- ``mark`` -- (default: None) if set to i, the vertex i is highlighted.
297
- ``save_pos`` -- (default:False) if True, the positions of the vertices are saved.
298
299
EXAMPLES::
300
301
sage: S = ClusterSeed(['F',4,[1,2]])
302
sage: S.save_image(os.path.join(SAGE_TMP, 'sage.png'))
303
"""
304
graph_plot = self.plot( circular=circular, mark=mark, save_pos=save_pos)
305
graph_plot.save( filename=filename )
306
307
def b_matrix(self):
308
r"""
309
Returns the `B` *-matrix* of ``self``.
310
311
EXAMPLES::
312
313
sage: ClusterSeed(['A',4]).b_matrix()
314
[ 0 1 0 0]
315
[-1 0 -1 0]
316
[ 0 1 0 1]
317
[ 0 0 -1 0]
318
319
sage: ClusterSeed(['B',4]).b_matrix()
320
[ 0 1 0 0]
321
[-1 0 -1 0]
322
[ 0 1 0 1]
323
[ 0 0 -2 0]
324
325
sage: ClusterSeed(['D',4]).b_matrix()
326
[ 0 1 0 0]
327
[-1 0 -1 -1]
328
[ 0 1 0 0]
329
[ 0 1 0 0]
330
331
sage: ClusterSeed(QuiverMutationType([['A',2],['B',2]])).b_matrix()
332
[ 0 1 0 0]
333
[-1 0 0 0]
334
[ 0 0 0 1]
335
[ 0 0 -2 0]
336
"""
337
return copy( self._M )
338
339
def ground_field(self):
340
r"""
341
Returns the *ground field* of the cluster of ``self``.
342
343
EXAMPLES::
344
345
sage: S = ClusterSeed(['A',3])
346
sage: S.ground_field()
347
Fraction Field of Multivariate Polynomial Ring in x0, x1, x2 over Rational Field
348
"""
349
return self._R
350
351
def x(self,k):
352
r"""
353
Returns the `k` *-th initial cluster variable* for the associated cluster seed.
354
355
EXAMPLES::
356
357
sage: S = ClusterSeed(['A',3])
358
sage: S.mutate([2,1])
359
sage: S.x(0)
360
x0
361
362
sage: S.x(1)
363
x1
364
365
sage: S.x(2)
366
x2
367
"""
368
if k in range(self._n):
369
x = self._R.gens()[k]
370
return ClusterVariable( x.parent(), x.numerator(), x.denominator(), mutation_type=self._mutation_type, variable_type='cluster variable' )
371
else:
372
raise ValueError("The input is not in an index of a cluster variable.")
373
374
def y(self,k):
375
r"""
376
Returns the `k` *-th initial coefficient (frozen variable)* for the associated cluster seed.
377
378
EXAMPLES::
379
380
sage: S = ClusterSeed(['A',3]).principal_extension()
381
sage: S.mutate([2,1])
382
sage: S.y(0)
383
y0
384
385
sage: S.y(1)
386
y1
387
388
sage: S.y(2)
389
y2
390
"""
391
if k in range(self._m):
392
x = self._R.gens()[self._n+k]
393
return ClusterVariable( x.parent(), x.numerator(), x.denominator(), mutation_type=self._mutation_type, variable_type='frozen variable' )
394
else:
395
raise ValueError("The input is not in an index of a frozen variable.")
396
397
def n(self):
398
r"""
399
Returns the number of *exchangeable variables* of ``self``.
400
401
EXAMPLES::
402
403
sage: S = ClusterSeed(['A',3])
404
sage: S.n()
405
3
406
"""
407
return self._n
408
409
def m(self):
410
r"""
411
Returns the number of *frozen variables* of ``self``.
412
413
EXAMPLES::
414
415
sage: S = ClusterSeed(['A',3])
416
sage: S.n()
417
3
418
419
sage: S.m()
420
0
421
422
sage: S = S.principal_extension()
423
sage: S.m()
424
3
425
"""
426
return self._m
427
428
def cluster_variable(self,k):
429
r"""
430
Returns the `k`-th *cluster variable* of ``self``.
431
432
EXAMPLES::
433
434
sage: S = ClusterSeed(['A',3])
435
sage: S.mutate([1,2])
436
437
sage: [S.cluster_variable(k) for k in range(3)]
438
[x0, (x0*x2 + 1)/x1, (x0*x2 + x1 + 1)/(x1*x2)]
439
"""
440
if k not in range(self._n):
441
raise ValueError("The cluster seed does not have a cluster variable of index %s."%k)
442
f = self._cluster[k]
443
return ClusterVariable( f.parent(), f.numerator(), f.denominator(), mutation_type=self._mutation_type, variable_type='cluster variable' )
444
445
def cluster(self):
446
r"""
447
Returns the *cluster* of ``self``.
448
449
EXAMPLES::
450
451
sage: S = ClusterSeed(['A',3])
452
sage: S.cluster()
453
[x0, x1, x2]
454
455
sage: S.mutate(1)
456
sage: S.cluster()
457
[x0, (x0*x2 + 1)/x1, x2]
458
459
sage: S.mutate(2)
460
sage: S.cluster()
461
[x0, (x0*x2 + 1)/x1, (x0*x2 + x1 + 1)/(x1*x2)]
462
463
sage: S.mutate([2,1])
464
sage: S.cluster()
465
[x0, x1, x2]
466
"""
467
return [ self.cluster_variable(k) for k in range(self._n) ]
468
469
def f_polynomial(self,k,ignore_coefficients=False):
470
r"""
471
Returns the ``k``-th *F-polynomial* of ``self``. It is obtained from the
472
``k``-th cluster variable by setting all `x_i` to `1`.
473
474
Requires principal coefficients, initialized by using principal_extension(),
475
or the user can set 'ignore_coefficients=True' to bypass this restriction.
476
477
Warning: this method assumes the sign-coherence conjecture and that the
478
input seed is sign-coherent (has an exchange matrix with columns of like signs).
479
Otherwise, computational errors might arise.
480
481
EXAMPLES::
482
483
sage: S = ClusterSeed(['A',3]).principal_extension()
484
sage: S.mutate([2,1,2])
485
sage: [S.f_polynomial(k) for k in range(3)]
486
[1, y1*y2 + y2 + 1, y1 + 1]
487
488
sage: S = ClusterSeed(Matrix([[0,1],[-1,0],[1,0],[-1,1]])); S
489
A seed for a cluster algebra of rank 2 with 2 frozen variables
490
sage: T = ClusterSeed(Matrix([[0,1],[-1,0]])).principal_extension(); T
491
A seed for a cluster algebra of rank 2 with principal coefficients
492
sage: S.mutate(0)
493
sage: T.mutate(0)
494
sage: S.f_polynomials()
495
Traceback (most recent call last):
496
...
497
ValueError: No principal coefficients initialized. Use principal_extension, or ignore_coefficients to ignore this.
498
sage: S.f_polynomials(ignore_coefficients=True)
499
[y0 + y1, 1]
500
sage: T.f_polynomials()
501
[y0 + 1, 1]
502
"""
503
if not (ignore_coefficients or self._is_principal):
504
raise ValueError("No principal coefficients initialized. Use principal_extension, or ignore_coefficients to ignore this.")
505
506
if k not in range(self._n):
507
raise ValueError("The cluster seed does not have a cluster variable of index %s."%k)
508
eval_dict = dict( [ ( self.x(i), 1 ) for i in range(self._n) ] )
509
return self.cluster_variable(k).subs(eval_dict)
510
511
def f_polynomials(self,ignore_coefficients=False):
512
r"""
513
Returns all *F-polynomials* of ``self``. These are obtained from the
514
cluster variables by setting all `x_i`'s to `1`.
515
516
Requires principal coefficients, initialized by using principal_extension(),
517
or the user can set 'ignore_coefficients=True' to bypass this restriction.
518
519
Warning: this method assumes the sign-coherence conjecture and that the
520
input seed is sign-coherent (has an exchange matrix with columns of like signs).
521
Otherwise, computational errors might arise.
522
523
EXAMPLES::
524
525
sage: S = ClusterSeed(['A',3]).principal_extension()
526
sage: S.mutate([2,1,2])
527
sage: S.f_polynomials()
528
[1, y1*y2 + y2 + 1, y1 + 1]
529
"""
530
if not (ignore_coefficients or self._is_principal):
531
raise ValueError("No principal coefficients initialized. Use principal_extension, or ignore_coefficients to ignore this.")
532
533
return [ self.f_polynomial(k,ignore_coefficients=ignore_coefficients) for k in range(self._n) ]
534
535
def g_vector(self,k,ignore_coefficients=False):
536
r"""
537
Returns the ``k``-th *g-vector* of ``self``. This is the degree vector
538
of the ``k``-th cluster variable after setting all `y_i`'s to `0`.
539
540
Requires principal coefficients, initialized by using principal_extension(),
541
or the user can set 'ignore_coefficients=True' to bypass this restriction.
542
543
Warning: this method assumes the sign-coherence conjecture and that the
544
input seed is sign-coherent (has an exchange matrix with columns of like signs).
545
Otherwise, computational errors might arise.
546
547
EXAMPLES::
548
549
sage: S = ClusterSeed(['A',3]).principal_extension()
550
sage: S.mutate([2,1,2])
551
sage: [ S.g_vector(k) for k in range(3) ]
552
[(1, 0, 0), (0, 0, -1), (0, -1, 0)]
553
"""
554
if not (ignore_coefficients or self._is_principal):
555
raise ValueError("No principal coefficients initialized. Use principal_extension, or ignore_coefficients to ignore this.")
556
if k not in range(self._n):
557
raise ValueError("The cluster seed does not have a cluster variable of index %s."%k)
558
f = self.cluster_variable(k)
559
eval_dict = dict( [ ( self.y(i), 0 ) for i in range(self._m) ] )
560
f0 = f.subs(eval_dict)
561
d1 = f0.numerator().degrees()
562
d2 = f0.denominator().degrees()
563
return tuple( d1[i] - d2[i] for i in range(self._n) )
564
565
def g_matrix(self,ignore_coefficients=False):
566
r"""
567
Returns the matrix of all *g-vectors* of ``self``. This are the degree vectors
568
of the cluster variables after setting all `y_i`'s to `0`.
569
570
Requires principal coefficients, initialized by using principal_extension(),
571
or the user can set 'ignore_coefficients=True' to bypass this restriction.
572
573
Warning: this method assumes the sign-coherence conjecture and that the
574
input seed is sign-coherent (has an exchange matrix with columns of like signs).
575
Otherwise, computational errors might arise.
576
577
EXAMPLES::
578
579
sage: S = ClusterSeed(['A',3]).principal_extension()
580
sage: S.mutate([2,1,2])
581
sage: S.g_matrix()
582
[ 1 0 0]
583
[ 0 0 -1]
584
[ 0 -1 0]
585
586
sage: S = ClusterSeed(['A',3])
587
sage: S2 = S.principal_extension()
588
sage: S.mutate([0,1])
589
sage: S2.mutate([0,1])
590
sage: S.g_matrix()
591
Traceback (most recent call last):
592
...
593
ValueError: No principal coefficients initialized. Use
594
principal_extension, or ignore_coefficients to ignore this.
595
sage: S.g_matrix(ignore_coefficients=True)
596
[-1 0 0]
597
[ 1 0 0]
598
[ 0 1 1]
599
sage: S2.g_matrix()
600
[-1 -1 0]
601
[ 1 0 0]
602
[ 0 0 1]
603
"""
604
from sage.matrix.all import matrix
605
if not (ignore_coefficients or self._is_principal):
606
raise ValueError("No principal coefficients initialized. Use principal_extension, or ignore_coefficients to ignore this.")
607
return matrix( [ self.g_vector(k,ignore_coefficients=ignore_coefficients) for k in range(self._n) ] ).transpose()
608
609
def c_vector(self,k,ignore_coefficients=False):
610
r"""
611
Returns the ``k``-th *c-vector* of ``self``. It is obtained as the
612
``k``-th column vector of the bottom part of the ``B``-matrix of ``self``.
613
614
Requires principal coefficients, initialized by using principal_extension(),
615
or the user can set 'ignore_coefficients=True' to bypass this restriction.
616
617
Warning: this method assumes the sign-coherence conjecture and that the
618
input seed is sign-coherent (has an exchange matrix with columns of like signs).
619
Otherwise, computational errors might arise.
620
621
EXAMPLES::
622
623
sage: S = ClusterSeed(['A',3]).principal_extension()
624
sage: S.mutate([2,1,2])
625
sage: [ S.c_vector(k) for k in range(3) ]
626
[(1, 0, 0), (0, 0, -1), (0, -1, 0)]
627
628
sage: S = ClusterSeed(Matrix([[0,1],[-1,0],[1,0],[-1,1]])); S
629
A seed for a cluster algebra of rank 2 with 2 frozen variables
630
sage: S.c_vector(0)
631
Traceback (most recent call last):
632
...
633
ValueError: No principal coefficients initialized. Use principal_extension, or ignore_coefficients to ignore this.
634
sage: S.c_vector(0,ignore_coefficients=True)
635
(1, -1)
636
"""
637
if k not in range(self._n):
638
raise ValueError("The cluster seed does not have a c-vector of index %s."%k)
639
if not (ignore_coefficients or self._is_principal):
640
raise ValueError("No principal coefficients initialized. Use principal_extension, or ignore_coefficients to ignore this.")
641
return tuple( self._M[i,k] for i in range(self._n,self._n+self._m) )
642
643
def c_matrix(self,ignore_coefficients=False):
644
r"""
645
Returns all *c-vectors* of ``self``.
646
647
Requires principal coefficients, initialized by using principal_extension(),
648
or the user can set 'ignore_coefficients=True' to bypass this restriction.
649
650
Warning: this method assumes the sign-coherence conjecture and that the
651
input seed is sign-coherent (has an exchange matrix with columns of like signs).
652
Otherwise, computational errors might arise.
653
654
EXAMPLES::
655
656
sage: S = ClusterSeed(['A',3]).principal_extension()
657
sage: S.mutate([2,1,2])
658
sage: S.c_matrix()
659
[ 1 0 0]
660
[ 0 0 -1]
661
[ 0 -1 0]
662
"""
663
if not (ignore_coefficients or self._is_principal):
664
raise ValueError("No principal coefficients initialized. Use principal_extension, or ignore_coefficients to ignore this.")
665
666
return self._M.submatrix(self._n,0)
667
668
def coefficient(self,k):
669
r"""
670
Returns the *coefficient* of ``self`` at index ``k``.
671
672
EXAMPLES::
673
674
sage: S = ClusterSeed(['A',3]).principal_extension()
675
sage: S.mutate([2,1,2])
676
sage: [ S.coefficient(k) for k in range(3) ]
677
[y0, 1/y2, 1/y1]
678
"""
679
from sage.misc.all import prod
680
if k not in range(self._n):
681
raise ValueError("The cluster seed does not have a coefficient of index %s."%k)
682
if self._m == 0:
683
return self.x(0)**0
684
#### Note: this special case m = 0 no longer needed except if we want type(answer) to be a cluster variable rather than an integer.
685
else:
686
exp = self.c_vector(k,ignore_coefficients=True)
687
return prod( self.y(i)**exp[i] for i in xrange(self._m) )
688
689
def coefficients(self):
690
r"""
691
Returns all *coefficients* of ``self``.
692
693
EXAMPLES::
694
695
sage: S = ClusterSeed(['A',3]).principal_extension()
696
sage: S.mutate([2,1,2])
697
sage: S.coefficients()
698
[y0, 1/y2, 1/y1]
699
"""
700
return [ self.coefficient(k) for k in range(self._n) ]
701
702
def quiver(self):
703
r"""
704
Returns the *quiver* associated to ``self``.
705
706
EXAMPLES::
707
708
sage: S = ClusterSeed(['A',3])
709
sage: S.quiver()
710
Quiver on 3 vertices of type ['A', 3]
711
"""
712
from sage.combinat.cluster_algebra_quiver.quiver import ClusterQuiver
713
if self._quiver is None:
714
self._quiver = ClusterQuiver( self._M )
715
return self._quiver
716
717
def is_acyclic(self):
718
r"""
719
Returns True iff self is acyclic (i.e., if the underlying quiver is acyclic).
720
721
EXAMPLES::
722
723
sage: ClusterSeed(['A',4]).is_acyclic()
724
True
725
726
sage: ClusterSeed(['A',[2,1],1]).is_acyclic()
727
True
728
729
sage: ClusterSeed([[0,1],[1,2],[2,0]]).is_acyclic()
730
False
731
"""
732
return self.quiver()._digraph.is_directed_acyclic()
733
734
def is_bipartite(self,return_bipartition=False):
735
r"""
736
Returns True iff self is bipartite (i.e., if the underlying quiver is bipartite).
737
738
INPUT:
739
740
- return_bipartition -- (default:False) if True, the bipartition is returned in the case of ``self`` being bipartite.
741
742
EXAMPLES::
743
744
sage: ClusterSeed(['A',[3,3],1]).is_bipartite()
745
True
746
747
sage: ClusterSeed(['A',[4,3],1]).is_bipartite()
748
False
749
"""
750
return self.quiver().is_bipartite(return_bipartition=return_bipartition)
751
752
def mutate(self, sequence, inplace=True):
753
r"""
754
Mutates ``self`` at a vertex or a sequence of vertices.
755
756
INPUT:
757
758
- ``sequence`` -- a vertex of self or an iterator of vertices of self.
759
- ``inplace`` -- (default: True) if False, the result is returned, otherwise ``self`` is modified.
760
761
EXAMPLES::
762
763
sage: S = ClusterSeed(['A',4]); S.b_matrix()
764
[ 0 1 0 0]
765
[-1 0 -1 0]
766
[ 0 1 0 1]
767
[ 0 0 -1 0]
768
769
sage: S.mutate(0); S.b_matrix()
770
[ 0 -1 0 0]
771
[ 1 0 -1 0]
772
[ 0 1 0 1]
773
[ 0 0 -1 0]
774
775
sage: T = S.mutate(0, inplace=False); T
776
A seed for a cluster algebra of rank 4 of type ['A', 4]
777
778
sage: S.mutate(0)
779
sage: S == T
780
True
781
782
sage: S.mutate([0,1,0])
783
sage: S.b_matrix()
784
[ 0 -1 1 0]
785
[ 1 0 0 0]
786
[-1 0 0 1]
787
[ 0 0 -1 0]
788
789
sage: S = ClusterSeed(QuiverMutationType([['A',1],['A',3]]))
790
sage: S.b_matrix()
791
[ 0 0 0 0]
792
[ 0 0 1 0]
793
[ 0 -1 0 -1]
794
[ 0 0 1 0]
795
796
sage: T = S.mutate(0,inplace=False)
797
sage: S == T
798
False
799
"""
800
if inplace:
801
seed = self
802
else:
803
seed = ClusterSeed( self )
804
805
n, m = seed._n, seed._m
806
V = range(n)
807
808
if sequence in V:
809
seq = [sequence]
810
else:
811
seq = sequence
812
if type( seq ) is tuple:
813
seq = list( seq )
814
if not type( seq ) is list:
815
raise ValueError('The quiver can only be mutated at a vertex or at a sequence of vertices')
816
if not type(inplace) is bool:
817
raise ValueError('The second parameter must be boolean. To mutate at a sequence of length 2, input it as a list.')
818
if any( v not in V for v in seq ):
819
v = filter( lambda v: v not in V, seq )[0]
820
raise ValueError('The quiver cannot be mutated at the vertex ' + str( v ))
821
822
for k in seq:
823
M = seed._M
824
cluster = seed._cluster
825
mon_p = seed._R(1)
826
mon_n = seed._R(1)
827
828
for j in range(n+m):
829
if M[j,k] > 0:
830
mon_p = mon_p*cluster[j]**M[j,k]
831
elif M[j,k] < 0:
832
mon_n = mon_n*cluster[j]**(-M[j,k])
833
834
cluster[k] = (mon_p+mon_n)*cluster[k]**(-1)
835
seed._M.mutate(k)
836
#seed._M = _matrix_mutate( seed._M, k )
837
838
seed._quiver = None
839
if not inplace:
840
return seed
841
842
def mutation_sequence(self, sequence, show_sequence=False, fig_size=1.2,return_output='seed'):
843
r"""
844
Returns the seeds obtained by mutating ``self`` at all vertices in ``sequence``.
845
846
INPUT:
847
848
- ``sequence`` -- an iterable of vertices of self.
849
- ``show_sequence`` -- (default: False) if True, a png containing the associated quivers is shown.
850
- ``fig_size`` -- (default: 1.2) factor by which the size of the plot is multiplied.
851
- ``return_output`` -- (default: 'seed') determines what output is to be returned::
852
853
* if 'seed', outputs all the cluster seeds obtained by the ``sequence`` of mutations.
854
* if 'matrix', outputs a list of exchange matrices.
855
* if 'var', outputs a list of new cluster variables obtained at each step.
856
857
EXAMPLES::
858
859
sage: S = ClusterSeed(['A',2])
860
sage: for T in S.mutation_sequence([0,1,0]):
861
... print T.b_matrix()
862
[ 0 -1]
863
[ 1 0]
864
[ 0 1]
865
[-1 0]
866
[ 0 -1]
867
[ 1 0]
868
869
sage: S=ClusterSeed(['A',2])
870
sage: S.mutation_sequence([0,1,0,1],return_output='var')
871
[(x1 + 1)/x0, (x0 + x1 + 1)/(x0*x1), (x0 + 1)/x1, x0]
872
"""
873
seed = ClusterSeed( self )
874
875
new_clust_var = []
876
seed_sequence = []
877
878
for v in sequence:
879
seed = seed.mutate(v,inplace=False)
880
new_clust_var.append( seed._cluster[v])
881
seed_sequence.append( seed )
882
883
if show_sequence:
884
self.quiver().mutation_sequence(sequence=sequence, show_sequence=True, fig_size=fig_size )
885
886
if return_output=='seed':
887
return seed_sequence
888
elif return_output=='matrix':
889
return [ seed._M for seed in seed_sequence ]
890
elif return_output=='var':
891
return new_clust_var
892
else:
893
raise ValueError('The parameter `return_output` can only be `seed`, `matrix`, or `var`.')
894
895
def exchangeable_part(self):
896
r"""
897
Returns the restriction to the principal part (i.e. the exchangeable variables) of ``self``.
898
899
EXAMPLES::
900
901
sage: S = ClusterSeed(['A',4])
902
sage: T = ClusterSeed( S.quiver().digraph().edges(), frozen=1 )
903
sage: T.quiver().digraph().edges()
904
[(0, 1, (1, -1)), (2, 1, (1, -1)), (2, 3, (1, -1))]
905
906
sage: T.exchangeable_part().quiver().digraph().edges()
907
[(0, 1, (1, -1)), (2, 1, (1, -1))]
908
909
sage: S2 = S.principal_extension()
910
sage: S3 = S2.principal_extension(ignore_coefficients=True)
911
sage: S2.exchangeable_part() == S3.exchangeable_part()
912
True
913
"""
914
from sage.combinat.cluster_algebra_quiver.mutation_class import _principal_part
915
eval_dict = dict( [ ( self.y(i), 1 ) for i in xrange(self._m) ] )
916
seed = ClusterSeed( _principal_part( self._M ) )
917
seed._cluster = [ self._cluster[k].subs(eval_dict) for k in xrange(self._n) ]
918
seed._mutation_type = self._mutation_type
919
return seed
920
921
def principal_extension(self,ignore_coefficients=False):
922
r"""
923
Returns the principal extension of self, yielding a 2n-by-n matrix. Raises an error if the input seed has a non-square exchange matrix,
924
unless 'ignore_coefficients=True' is set. In this case, the method instead adds n frozen variables to any previously frozen variables.
925
I.e., the seed obtained by adding a frozen variable to every exchangeable variable of ``self``.
926
927
EXAMPLES::
928
929
sage: S = ClusterSeed([[0,1],[1,2],[2,3],[2,4]]); S
930
A seed for a cluster algebra of rank 5
931
932
sage: T = S.principal_extension(); T
933
A seed for a cluster algebra of rank 5 with principal coefficients
934
935
sage: T.b_matrix()
936
[ 0 1 0 0 0]
937
[-1 0 1 0 0]
938
[ 0 -1 0 1 1]
939
[ 0 0 -1 0 0]
940
[ 0 0 -1 0 0]
941
[ 1 0 0 0 0]
942
[ 0 1 0 0 0]
943
[ 0 0 1 0 0]
944
[ 0 0 0 1 0]
945
[ 0 0 0 0 1]
946
947
sage: T2 = T.principal_extension()
948
Traceback (most recent call last):
949
...
950
ValueError: The b-matrix is not square. Use ignore_coefficients to ignore this.
951
952
sage: T2 = T.principal_extension(ignore_coefficients=True); T2.b_matrix()
953
[ 0 1 0 0 0]
954
[-1 0 1 0 0]
955
[ 0 -1 0 1 1]
956
[ 0 0 -1 0 0]
957
[ 0 0 -1 0 0]
958
[ 1 0 0 0 0]
959
[ 0 1 0 0 0]
960
[ 0 0 1 0 0]
961
[ 0 0 0 1 0]
962
[ 0 0 0 0 1]
963
[ 1 0 0 0 0]
964
[ 0 1 0 0 0]
965
[ 0 0 1 0 0]
966
[ 0 0 0 1 0]
967
[ 0 0 0 0 1]
968
"""
969
from sage.matrix.all import identity_matrix
970
if not ignore_coefficients and self._m != 0:
971
raise ValueError("The b-matrix is not square. Use ignore_coefficients to ignore this.")
972
M = self._M.stack(identity_matrix(self._n))
973
is_principal = (self._m == 0)
974
seed = ClusterSeed( M, is_principal=is_principal )
975
seed._mutation_type = self._mutation_type
976
return seed
977
978
def reorient( self, data ):
979
r"""
980
Reorients ``self`` with respect to the given total order,
981
or with respect to an iterator of ordered pairs.
982
983
WARNING:
984
985
- This operation might change the mutation type of ``self``.
986
- Ignores ordered pairs `(i,j)` for which neither `(i,j)` nor `(j,i)` is an edge of ``self``.
987
988
INPUT:
989
990
- ``data`` -- an iterator defining a total order on ``self.vertices()``, or an iterator of ordered pairs in ``self`` defining the new orientation of these edges.
991
992
EXAMPLES::
993
994
sage: S = ClusterSeed(['A',[2,3],1])
995
sage: S.mutation_type()
996
['A', [2, 3], 1]
997
998
sage: S.reorient([(0,1),(2,3)])
999
sage: S.mutation_type()
1000
['D', 5]
1001
1002
sage: S.reorient([(1,0),(2,3)])
1003
sage: S.mutation_type()
1004
['A', [1, 4], 1]
1005
1006
sage: S.reorient([0,1,2,3,4])
1007
sage: S.mutation_type()
1008
['A', [1, 4], 1]
1009
"""
1010
if not self._quiver:
1011
self.quiver()
1012
self._quiver.reorient( data )
1013
self._M = self._quiver._M
1014
self.reset_cluster()
1015
self._mutation_type = None
1016
1017
def set_cluster( self, cluster ):
1018
r"""
1019
Sets the cluster for ``self`` to ``cluster``.
1020
1021
INPUT:
1022
1023
- ``cluster`` -- an iterable defining a cluster for ``self``.
1024
1025
EXAMPLES::
1026
1027
sage: S = ClusterSeed(['A',3])
1028
sage: cluster = S.cluster()
1029
sage: S.mutate([1,2,1])
1030
sage: S.cluster()
1031
[x0, (x1 + 1)/x2, (x0*x2 + x1 + 1)/(x1*x2)]
1032
1033
sage: S.set_cluster(cluster)
1034
sage: S.cluster()
1035
[x0, x1, x2]
1036
"""
1037
if not len(cluster) == self._n+self._m:
1038
raise ValueError('The number of given cluster variables is wrong')
1039
if any(c not in self._R for c in cluster):
1040
raise ValueError('The cluster variables are not all contained in %s'%self._R)
1041
self._cluster = [ self._R(x) for x in cluster ]
1042
self._is_principal = None
1043
1044
def reset_cluster( self ):
1045
r"""
1046
Resets the cluster of ``self`` to the initial cluster.
1047
1048
EXAMPLES::
1049
1050
sage: S = ClusterSeed(['A',3])
1051
sage: S.mutate([1,2,1])
1052
sage: S.cluster()
1053
[x0, (x1 + 1)/x2, (x0*x2 + x1 + 1)/(x1*x2)]
1054
1055
sage: S.reset_cluster()
1056
sage: S.cluster()
1057
[x0, x1, x2]
1058
1059
sage: T = S.principal_extension()
1060
sage: T.cluster()
1061
[x0, x1, x2]
1062
sage: T.mutate([1,2,1])
1063
sage: T.cluster()
1064
[x0, (x1*y2 + x0)/x2, (x1*y1*y2 + x0*y1 + x2)/(x1*x2)]
1065
1066
sage: T.reset_cluster()
1067
sage: T.cluster()
1068
[x0, x1, x2]
1069
"""
1070
self.set_cluster(self._R.gens())
1071
1072
def reset_coefficients( self ):
1073
r"""
1074
Resets the coefficients of ``self`` to the frozen variables but keeps the current cluster.
1075
Raises an error if the number of frozen variables is different than the number of exchangeable variables.
1076
1077
EXAMPLES::
1078
1079
sage: S = ClusterSeed(['A',3]).principal_extension()
1080
sage: S.b_matrix()
1081
[ 0 1 0]
1082
[-1 0 -1]
1083
[ 0 1 0]
1084
[ 1 0 0]
1085
[ 0 1 0]
1086
[ 0 0 1]
1087
sage: S.mutate([1,2,1])
1088
sage: S.b_matrix()
1089
[ 0 1 -1]
1090
[-1 0 1]
1091
[ 1 -1 0]
1092
[ 1 0 0]
1093
[ 0 1 -1]
1094
[ 0 0 -1]
1095
sage: S.reset_coefficients()
1096
sage: S.b_matrix()
1097
[ 0 1 -1]
1098
[-1 0 1]
1099
[ 1 -1 0]
1100
[ 1 0 0]
1101
[ 0 1 0]
1102
[ 0 0 1]
1103
"""
1104
n,m = self._n, self._m
1105
if not n == m:
1106
raise ValueError("The numbers of cluster variables and of frozen variables do not coincide.")
1107
for i in xrange(m):
1108
for j in xrange(n):
1109
if i == j:
1110
self._M[i+n,j] = 1
1111
else:
1112
self._M[i+n,j] = 0
1113
self._quiver = None
1114
self._is_principal = None
1115
1116
def mutation_class_iter( self, depth=infinity, show_depth=False, return_paths=False, up_to_equivalence=True, only_sink_source=False ):
1117
r"""
1118
Returns an iterator for the mutation class of ``self`` with respect to certain constrains.
1119
1120
INPUT:
1121
1122
- ``depth`` -- (default: infinity) integer or infinity, only seeds with distance at most ``depth`` from ``self`` are returned.
1123
- ``show_depth`` -- (default: False) if True, the current depth of the mutation is shown while computing.
1124
- ``return_paths`` -- (default: False) if True, a shortest path of mutations from ``self`` to the given quiver is returned as well.
1125
- ``up_to_equivalence`` -- (default: True) if True, only one seed up to simultaneous permutation of rows and columns of the exchange matrix is recorded.
1126
- ``sink_source`` -- (default: False) if True, only mutations at sinks and sources are applied.
1127
1128
EXAMPLES:
1129
1130
A standard finite type example::
1131
1132
sage: S = ClusterSeed(['A',3])
1133
sage: it = S.mutation_class_iter()
1134
sage: for T in it: print T
1135
A seed for a cluster algebra of rank 3 of type ['A', 3]
1136
A seed for a cluster algebra of rank 3 of type ['A', 3]
1137
A seed for a cluster algebra of rank 3 of type ['A', 3]
1138
A seed for a cluster algebra of rank 3 of type ['A', 3]
1139
A seed for a cluster algebra of rank 3 of type ['A', 3]
1140
A seed for a cluster algebra of rank 3 of type ['A', 3]
1141
A seed for a cluster algebra of rank 3 of type ['A', 3]
1142
A seed for a cluster algebra of rank 3 of type ['A', 3]
1143
A seed for a cluster algebra of rank 3 of type ['A', 3]
1144
A seed for a cluster algebra of rank 3 of type ['A', 3]
1145
A seed for a cluster algebra of rank 3 of type ['A', 3]
1146
A seed for a cluster algebra of rank 3 of type ['A', 3]
1147
A seed for a cluster algebra of rank 3 of type ['A', 3]
1148
A seed for a cluster algebra of rank 3 of type ['A', 3]
1149
1150
A finite type example with given depth::
1151
1152
sage: it = S.mutation_class_iter(depth=1)
1153
sage: for T in it: print T
1154
A seed for a cluster algebra of rank 3 of type ['A', 3]
1155
A seed for a cluster algebra of rank 3 of type ['A', 3]
1156
A seed for a cluster algebra of rank 3 of type ['A', 3]
1157
A seed for a cluster algebra of rank 3 of type ['A', 3]
1158
1159
A finite type example where the depth is shown while computing::
1160
1161
sage: it = S.mutation_class_iter(show_depth=True)
1162
sage: for T in it: pass
1163
Depth: 0 found: 1 Time: ... s
1164
Depth: 1 found: 4 Time: ... s
1165
Depth: 2 found: 9 Time: ... s
1166
Depth: 3 found: 13 Time: ... s
1167
Depth: 4 found: 14 Time: ... s
1168
1169
A finite type example with shortest paths returned::
1170
1171
sage: it = S.mutation_class_iter(return_paths=True)
1172
sage: for T in it: print T
1173
(A seed for a cluster algebra of rank 3 of type ['A', 3], [])
1174
(A seed for a cluster algebra of rank 3 of type ['A', 3], [2])
1175
(A seed for a cluster algebra of rank 3 of type ['A', 3], [1])
1176
(A seed for a cluster algebra of rank 3 of type ['A', 3], [0])
1177
(A seed for a cluster algebra of rank 3 of type ['A', 3], [2, 1])
1178
(A seed for a cluster algebra of rank 3 of type ['A', 3], [0, 2])
1179
(A seed for a cluster algebra of rank 3 of type ['A', 3], [0, 1])
1180
(A seed for a cluster algebra of rank 3 of type ['A', 3], [1, 2])
1181
(A seed for a cluster algebra of rank 3 of type ['A', 3], [1, 0])
1182
(A seed for a cluster algebra of rank 3 of type ['A', 3], [0, 2, 1])
1183
(A seed for a cluster algebra of rank 3 of type ['A', 3], [0, 1, 2])
1184
(A seed for a cluster algebra of rank 3 of type ['A', 3], [2, 1, 0])
1185
(A seed for a cluster algebra of rank 3 of type ['A', 3], [1, 0, 2])
1186
(A seed for a cluster algebra of rank 3 of type ['A', 3], [0, 1, 2, 0])
1187
1188
Finite type examples not considered up to equivalence::
1189
1190
sage: it = S.mutation_class_iter(up_to_equivalence=False)
1191
sage: len( [ T for T in it ] )
1192
84
1193
1194
sage: it = ClusterSeed(['A',2]).mutation_class_iter(return_paths=True,up_to_equivalence=False)
1195
sage: for T in it: print T
1196
(A seed for a cluster algebra of rank 2 of type ['A', 2], [])
1197
(A seed for a cluster algebra of rank 2 of type ['A', 2], [1])
1198
(A seed for a cluster algebra of rank 2 of type ['A', 2], [0])
1199
(A seed for a cluster algebra of rank 2 of type ['A', 2], [0, 1])
1200
(A seed for a cluster algebra of rank 2 of type ['A', 2], [1, 0])
1201
(A seed for a cluster algebra of rank 2 of type ['A', 2], [1, 0, 1])
1202
(A seed for a cluster algebra of rank 2 of type ['A', 2], [0, 1, 0])
1203
(A seed for a cluster algebra of rank 2 of type ['A', 2], [1, 0, 1, 0])
1204
(A seed for a cluster algebra of rank 2 of type ['A', 2], [0, 1, 0, 1])
1205
(A seed for a cluster algebra of rank 2 of type ['A', 2], [1, 0, 1, 0, 1])
1206
1207
Check that :trac:`14638` is fixed::
1208
1209
sage: S = ClusterSeed(['E',6])
1210
sage: MC = S.mutation_class(depth=7); len(MC)
1211
534
1212
1213
Infinite type examples::
1214
1215
sage: S = ClusterSeed(['A',[1,1],1])
1216
sage: it = S.mutation_class_iter()
1217
sage: it.next()
1218
A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1]
1219
sage: it.next()
1220
A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1]
1221
sage: it.next()
1222
A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1]
1223
sage: it.next()
1224
A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1]
1225
1226
sage: it = S.mutation_class_iter(depth=3, return_paths=True)
1227
sage: for T in it: print T
1228
(A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1], [])
1229
(A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1], [1])
1230
(A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1], [0])
1231
(A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1], [1, 0])
1232
(A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1], [0, 1])
1233
(A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1], [1, 0, 1])
1234
(A seed for a cluster algebra of rank 2 of type ['A', [1, 1], 1], [0, 1, 0])
1235
"""
1236
depth_counter = 0
1237
n = self._n
1238
timer = time.time()
1239
if return_paths:
1240
yield (self,[])
1241
else:
1242
yield self
1243
if up_to_equivalence:
1244
cl = Set( self._cluster )
1245
else:
1246
cl = tuple( self._cluster )
1247
clusters = {}
1248
clusters[ cl ] = [ self, range(n), [] ]
1249
gets_bigger = True
1250
if show_depth:
1251
timer2 = time.time()
1252
dc = str(depth_counter)
1253
dc += ' ' * (5-len(dc))
1254
nr = str(len(clusters))
1255
nr += ' ' * (10-len(nr))
1256
print "Depth: %s found: %s Time: %.2f s"%(dc,nr,timer2-timer)
1257
while gets_bigger and depth_counter < depth:
1258
gets_bigger = False
1259
keys = clusters.keys()
1260
for key in keys:
1261
sd = clusters[key]
1262
while sd[1]:
1263
i = sd[1].pop()
1264
if not only_sink_source or all( entry >= 0 for entry in sd[0]._M.row( i ) ) or all( entry <= 0 for entry in sd[0]._M.row( i ) ):
1265
sd2 = sd[0].mutate( i, inplace=False )
1266
if up_to_equivalence:
1267
cl2 = Set(sd2._cluster)
1268
else:
1269
cl2 = tuple(sd2._cluster)
1270
if cl2 in clusters:
1271
if not up_to_equivalence and i in clusters[cl2][1]:
1272
clusters[cl2][1].remove(i)
1273
else:
1274
gets_bigger = True
1275
if only_sink_source:
1276
orbits = range(n)
1277
else:
1278
orbits = [ index for index in xrange(n) if index > i or sd2._M[index,i] != 0 ]
1279
1280
clusters[ cl2 ] = [ sd2, orbits, clusters[key][2]+[i] ]
1281
if return_paths:
1282
yield (sd2,clusters[cl2][2])
1283
else:
1284
yield sd2
1285
depth_counter += 1
1286
if show_depth and gets_bigger:
1287
timer2 = time.time()
1288
dc = str(depth_counter)
1289
dc += ' ' * (5-len(dc))
1290
nr = str(len(clusters))
1291
nr += ' ' * (10-len(nr))
1292
print "Depth: %s found: %s Time: %.2f s"%(dc,nr,timer2-timer)
1293
1294
def mutation_class( self, depth=infinity, show_depth=False, return_paths=False, up_to_equivalence=True, only_sink_source=False ):
1295
r"""
1296
Returns the mutation class of ``self`` with respect to certain constraints.
1297
1298
INPUT:
1299
1300
- ``depth`` -- (default: infinity) integer, only seeds with distance at most depth from self are returned.
1301
- ``show_depth`` -- (default: False) if True, the actual depth of the mutation is shown.
1302
- ``return_paths`` -- (default: False) if True, a shortest path of mutation sequences from self to the given quiver is returned as well.
1303
- ``up_to_equivalence`` -- (default: True) if True, only seeds up to equivalence are considered.
1304
- ``sink_source`` -- (default: False) if True, only mutations at sinks and sources are applied.
1305
1306
EXAMPLES:
1307
1308
- for examples see :meth:`mutation_class_iter`
1309
1310
TESTS::
1311
1312
sage: A = ClusterSeed(['A',3]).mutation_class()
1313
"""
1314
if depth is infinity and not self.is_finite():
1315
raise ValueError('The mutation class can - for infinite types - only be computed up to a given depth')
1316
return list( S for S in self.mutation_class_iter( depth=depth, show_depth=show_depth, return_paths=return_paths, up_to_equivalence=up_to_equivalence, only_sink_source=only_sink_source ) )
1317
1318
def cluster_class_iter(self, depth=infinity, show_depth=False, up_to_equivalence=True):
1319
r"""
1320
Returns an iterator through all clusters in the mutation class of ``self``.
1321
1322
INPUT:
1323
1324
- ``depth`` -- (default: infinity) integer or infinity, only seeds with distance at most depth from self are returned
1325
- ``show_depth`` -- (default False) - if True, ignored if depth is set; returns the depth of the mutation class, i.e., the maximal distance from self of an element in the mutation class
1326
- ``up_to_equivalence`` -- (default: True) if True, only clusters up to equivalence are considered.
1327
1328
EXAMPLES:
1329
1330
A standard finite type example::
1331
1332
sage: S = ClusterSeed(['A',3])
1333
sage: it = S.cluster_class_iter()
1334
sage: for T in it: print T
1335
[x0, x1, x2]
1336
[x0, x1, (x1 + 1)/x2]
1337
[x0, (x0*x2 + 1)/x1, x2]
1338
[(x1 + 1)/x0, x1, x2]
1339
[x0, (x0*x2 + x1 + 1)/(x1*x2), (x1 + 1)/x2]
1340
[(x1 + 1)/x0, x1, (x1 + 1)/x2]
1341
[(x1 + 1)/x0, (x0*x2 + x1 + 1)/(x0*x1), x2]
1342
[x0, (x0*x2 + 1)/x1, (x0*x2 + x1 + 1)/(x1*x2)]
1343
[(x0*x2 + x1 + 1)/(x0*x1), (x0*x2 + 1)/x1, x2]
1344
[(x1 + 1)/x0, (x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2), (x1 + 1)/x2]
1345
[(x1 + 1)/x0, (x0*x2 + x1 + 1)/(x0*x1), (x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2)]
1346
[(x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2), (x0*x2 + x1 + 1)/(x1*x2), (x1 + 1)/x2]
1347
[(x0*x2 + x1 + 1)/(x0*x1), (x0*x2 + 1)/x1, (x0*x2 + x1 + 1)/(x1*x2)]
1348
[(x0*x2 + x1 + 1)/(x1*x2), (x0*x2 + x1 + 1)/(x0*x1), (x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2)]
1349
1350
A finite type example with given depth::
1351
1352
sage: it = S.cluster_class_iter(depth=1)
1353
sage: for T in it: print T
1354
[x0, x1, x2]
1355
[x0, x1, (x1 + 1)/x2]
1356
[x0, (x0*x2 + 1)/x1, x2]
1357
[(x1 + 1)/x0, x1, x2]
1358
1359
A finite type example where the depth is returned while computing::
1360
1361
sage: it = S.cluster_class_iter(show_depth=True)
1362
sage: for T in it: print T
1363
[x0, x1, x2]
1364
Depth: 0 found: 1 Time: ... s
1365
[x0, x1, (x1 + 1)/x2]
1366
[x0, (x0*x2 + 1)/x1, x2]
1367
[(x1 + 1)/x0, x1, x2]
1368
Depth: 1 found: 4 Time: ... s
1369
[x0, (x0*x2 + x1 + 1)/(x1*x2), (x1 + 1)/x2]
1370
[(x1 + 1)/x0, x1, (x1 + 1)/x2]
1371
[(x1 + 1)/x0, (x0*x2 + x1 + 1)/(x0*x1), x2]
1372
[x0, (x0*x2 + 1)/x1, (x0*x2 + x1 + 1)/(x1*x2)]
1373
[(x0*x2 + x1 + 1)/(x0*x1), (x0*x2 + 1)/x1, x2]
1374
Depth: 2 found: 9 Time: ... s
1375
[(x1 + 1)/x0, (x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2), (x1 + 1)/x2]
1376
[(x1 + 1)/x0, (x0*x2 + x1 + 1)/(x0*x1), (x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2)]
1377
[(x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2), (x0*x2 + x1 + 1)/(x1*x2), (x1 + 1)/x2]
1378
[(x0*x2 + x1 + 1)/(x0*x1), (x0*x2 + 1)/x1, (x0*x2 + x1 + 1)/(x1*x2)]
1379
Depth: 3 found: 13 Time: ... s
1380
[(x0*x2 + x1 + 1)/(x1*x2), (x0*x2 + x1 + 1)/(x0*x1), (x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2)]
1381
Depth: 4 found: 14 Time: ... s
1382
1383
Finite type examples not considered up to equivalence::
1384
1385
sage: it = S.cluster_class_iter(up_to_equivalence=False)
1386
sage: len( [ T for T in it ] )
1387
84
1388
1389
sage: it = ClusterSeed(['A',2]).cluster_class_iter(up_to_equivalence=False)
1390
sage: for T in it: print T
1391
[x0, x1]
1392
[x0, (x0 + 1)/x1]
1393
[(x1 + 1)/x0, x1]
1394
[(x1 + 1)/x0, (x0 + x1 + 1)/(x0*x1)]
1395
[(x0 + x1 + 1)/(x0*x1), (x0 + 1)/x1]
1396
[(x0 + x1 + 1)/(x0*x1), (x1 + 1)/x0]
1397
[(x0 + 1)/x1, (x0 + x1 + 1)/(x0*x1)]
1398
[x1, (x1 + 1)/x0]
1399
[(x0 + 1)/x1, x0]
1400
[x1, x0]
1401
1402
Infinite type examples::
1403
1404
sage: S = ClusterSeed(['A',[1,1],1])
1405
sage: it = S.cluster_class_iter()
1406
sage: it.next()
1407
[x0, x1]
1408
sage: it.next()
1409
[x0, (x0^2 + 1)/x1]
1410
sage: it.next()
1411
[(x1^2 + 1)/x0, x1]
1412
sage: it.next()
1413
[(x0^4 + 2*x0^2 + x1^2 + 1)/(x0*x1^2), (x0^2 + 1)/x1]
1414
sage: it.next()
1415
[(x1^2 + 1)/x0, (x1^4 + x0^2 + 2*x1^2 + 1)/(x0^2*x1)]
1416
1417
sage: it = S.cluster_class_iter(depth=3)
1418
sage: for T in it: print T
1419
[x0, x1]
1420
[x0, (x0^2 + 1)/x1]
1421
[(x1^2 + 1)/x0, x1]
1422
[(x0^4 + 2*x0^2 + x1^2 + 1)/(x0*x1^2), (x0^2 + 1)/x1]
1423
[(x1^2 + 1)/x0, (x1^4 + x0^2 + 2*x1^2 + 1)/(x0^2*x1)]
1424
[(x0^4 + 2*x0^2 + x1^2 + 1)/(x0*x1^2), (x0^6 + 3*x0^4 + 2*x0^2*x1^2 + x1^4 + 3*x0^2 + 2*x1^2 + 1)/(x0^2*x1^3)]
1425
[(x1^6 + x0^4 + 2*x0^2*x1^2 + 3*x1^4 + 2*x0^2 + 3*x1^2 + 1)/(x0^3*x1^2), (x1^4 + x0^2 + 2*x1^2 + 1)/(x0^2*x1)]
1426
"""
1427
mc_iter = self.mutation_class_iter( depth=depth, show_depth=show_depth, up_to_equivalence=up_to_equivalence )
1428
for c in mc_iter:
1429
yield c.cluster()
1430
1431
def cluster_class(self, depth=infinity, show_depth=False, up_to_equivalence=True):
1432
r"""
1433
Returns the cluster class of ``self`` with respect to certain constraints.
1434
1435
INPUT:
1436
1437
- ``depth`` -- (default: infinity) integer, only seeds with distance at most depth from self are returned
1438
- ``return_depth`` -- (default False) - if True, ignored if depth is set; returns the depth of the mutation class, i.e., the maximal distance from self of an element in the mutation class
1439
- ``up_to_equivalence`` -- (default: True) if True, only clusters up to equivalence are considered.
1440
1441
EXAMPLES:
1442
1443
- for examples see :meth:`cluster_class_iter`
1444
1445
TESTS::
1446
1447
sage: A = ClusterSeed(['A',3]).cluster_class()
1448
"""
1449
if depth is infinity and not self.is_finite():
1450
raise ValueError('The variable class can - for infinite types - only be computed up to a given depth')
1451
1452
return [ c for c in self.cluster_class_iter(depth=depth, show_depth=show_depth, up_to_equivalence=up_to_equivalence) ]
1453
1454
def b_matrix_class_iter(self, depth=infinity, up_to_equivalence=True):
1455
r"""
1456
Returns an iterator through all `B`-matrices in the mutation class of ``self``.
1457
1458
INPUT:
1459
1460
- ``depth`` -- (default:infinity) integer or infinity, only seeds with distance at most depth from self are returned
1461
- ``up_to_equivalence`` -- (default: True) if True, only 'B'-matrices up to equivalence are considered.
1462
1463
EXAMPLES:
1464
1465
A standard finite type example::
1466
1467
sage: S = ClusterSeed(['A',4])
1468
sage: it = S.b_matrix_class_iter()
1469
sage: for T in it: print T
1470
[ 0 0 0 1]
1471
[ 0 0 1 1]
1472
[ 0 -1 0 0]
1473
[-1 -1 0 0]
1474
[ 0 0 0 1]
1475
[ 0 0 1 0]
1476
[ 0 -1 0 1]
1477
[-1 0 -1 0]
1478
[ 0 0 1 1]
1479
[ 0 0 0 -1]
1480
[-1 0 0 0]
1481
[-1 1 0 0]
1482
[ 0 0 0 1]
1483
[ 0 0 -1 1]
1484
[ 0 1 0 -1]
1485
[-1 -1 1 0]
1486
[ 0 0 0 1]
1487
[ 0 0 -1 0]
1488
[ 0 1 0 -1]
1489
[-1 0 1 0]
1490
[ 0 0 0 -1]
1491
[ 0 0 -1 1]
1492
[ 0 1 0 -1]
1493
[ 1 -1 1 0]
1494
1495
A finite type example with given depth::
1496
1497
sage: it = S.b_matrix_class_iter(depth=1)
1498
sage: for T in it: print T
1499
[ 0 0 0 1]
1500
[ 0 0 1 1]
1501
[ 0 -1 0 0]
1502
[-1 -1 0 0]
1503
[ 0 0 0 1]
1504
[ 0 0 1 0]
1505
[ 0 -1 0 1]
1506
[-1 0 -1 0]
1507
[ 0 0 1 1]
1508
[ 0 0 0 -1]
1509
[-1 0 0 0]
1510
[-1 1 0 0]
1511
1512
Finite type example not considered up to equivalence::
1513
1514
sage: S = ClusterSeed(['A',3])
1515
sage: it = S.b_matrix_class_iter(up_to_equivalence=False)
1516
sage: for T in it: print T
1517
[ 0 1 0]
1518
[-1 0 -1]
1519
[ 0 1 0]
1520
[ 0 1 0]
1521
[-1 0 1]
1522
[ 0 -1 0]
1523
[ 0 -1 0]
1524
[ 1 0 1]
1525
[ 0 -1 0]
1526
[ 0 -1 0]
1527
[ 1 0 -1]
1528
[ 0 1 0]
1529
[ 0 -1 1]
1530
[ 1 0 -1]
1531
[-1 1 0]
1532
[ 0 1 -1]
1533
[-1 0 1]
1534
[ 1 -1 0]
1535
[ 0 0 1]
1536
[ 0 0 -1]
1537
[-1 1 0]
1538
[ 0 -1 1]
1539
[ 1 0 0]
1540
[-1 0 0]
1541
[ 0 0 -1]
1542
[ 0 0 1]
1543
[ 1 -1 0]
1544
[ 0 1 -1]
1545
[-1 0 0]
1546
[ 1 0 0]
1547
[ 0 1 1]
1548
[-1 0 0]
1549
[-1 0 0]
1550
[ 0 -1 -1]
1551
[ 1 0 0]
1552
[ 1 0 0]
1553
[ 0 0 -1]
1554
[ 0 0 -1]
1555
[ 1 1 0]
1556
[ 0 0 1]
1557
[ 0 0 1]
1558
[-1 -1 0]
1559
1560
Infinite (but finite mutation) type example::
1561
1562
sage: S = ClusterSeed(['A',[1,2],1])
1563
sage: it = S.b_matrix_class_iter()
1564
sage: for T in it: print T
1565
[ 0 1 1]
1566
[-1 0 1]
1567
[-1 -1 0]
1568
[ 0 -2 1]
1569
[ 2 0 -1]
1570
[-1 1 0]
1571
1572
Infinite mutation type example::
1573
1574
sage: S = ClusterSeed(['E',10])
1575
sage: it = S.b_matrix_class_iter(depth=3)
1576
sage: len ( [T for T in it] )
1577
266
1578
"""
1579
Q = self.quiver()
1580
for M in Q.mutation_class_iter( depth=depth, up_to_equivalence=up_to_equivalence, data_type='matrix' ):
1581
yield M
1582
1583
def b_matrix_class(self, depth=infinity, up_to_equivalence=True):
1584
r"""
1585
Returns all `B`-matrices in the mutation class of ``self``.
1586
1587
INPUT:
1588
1589
- ``depth`` -- (default:infinity) integer or infinity, only seeds with distance at most depth from self are returned
1590
- ``up_to_equivalence`` -- (default: True) if True, only 'B'-matrices up to equivalence are considered.
1591
1592
EXAMPLES:
1593
1594
- for examples see :meth:`b_matrix_class_iter`
1595
1596
TESTS::
1597
1598
sage: A = ClusterSeed(['A',3]).b_matrix_class()
1599
sage: A = ClusterSeed(['A',[2,1],1]).b_matrix_class()
1600
"""
1601
if depth is infinity and not self.is_mutation_finite():
1602
raise ValueError('The B-matrix class can - for infinite mutation types - only be computed up to a given depth')
1603
1604
return [ M for M in self.b_matrix_class_iter( depth=depth, up_to_equivalence=up_to_equivalence ) ]
1605
1606
def variable_class_iter(self, depth=infinity, ignore_bipartite_belt=False):
1607
r"""
1608
Returns an iterator for all cluster variables in the mutation class of ``self``.
1609
1610
INPUT:
1611
1612
- ``depth`` -- (default:infinity) integer, only seeds with distance at most depth from self are returned
1613
- ``ignore_bipartite_belt`` -- (default:False) if True, the algorithms does not use the bipartite belt
1614
1615
EXAMPLES:
1616
1617
A standard finite type example::
1618
1619
sage: S = ClusterSeed(['A',3])
1620
sage: it = S.variable_class_iter()
1621
sage: for T in it: print T
1622
x0
1623
x1
1624
x2
1625
(x1 + 1)/x0
1626
(x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2)
1627
(x1 + 1)/x2
1628
(x0*x2 + x1 + 1)/(x0*x1)
1629
(x0*x2 + 1)/x1
1630
(x0*x2 + x1 + 1)/(x1*x2)
1631
1632
Finite type examples with given depth::
1633
1634
sage: it = S.variable_class_iter(depth=1)
1635
sage: for T in it: print T
1636
Found a bipartite seed - restarting the depth counter at zero and constructing the variable class using its bipartite belt.
1637
x0
1638
x1
1639
x2
1640
(x1 + 1)/x0
1641
(x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2)
1642
(x1 + 1)/x2
1643
(x0*x2 + x1 + 1)/(x0*x1)
1644
(x0*x2 + 1)/x1
1645
(x0*x2 + x1 + 1)/(x1*x2)
1646
1647
Note that the notion of *depth* depends on whether a bipartite seed is found or not, or if it is manually ignored::
1648
1649
sage: it = S.variable_class_iter(depth=1,ignore_bipartite_belt=True)
1650
sage: for T in it: print T
1651
x0
1652
x1
1653
x2
1654
(x1 + 1)/x2
1655
(x0*x2 + 1)/x1
1656
(x1 + 1)/x0
1657
1658
sage: S.mutate([0,1])
1659
sage: it2 = S.variable_class_iter(depth=1)
1660
sage: for T in it2: print T
1661
(x1 + 1)/x0
1662
(x0*x2 + x1 + 1)/(x0*x1)
1663
x2
1664
(x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2)
1665
x1
1666
(x0*x2 + 1)/x1
1667
1668
Infinite type examples::
1669
1670
sage: S = ClusterSeed(['A',[1,1],1])
1671
sage: it = S.variable_class_iter(depth=2)
1672
sage: for T in it: print T
1673
Found a bipartite seed - restarting the depth counter at zero and constructing the variable class using its bipartite belt.
1674
x0
1675
x1
1676
(x1^2 + 1)/x0
1677
(x1^4 + x0^2 + 2*x1^2 + 1)/(x0^2*x1)
1678
(x0^4 + 2*x0^2 + x1^2 + 1)/(x0*x1^2)
1679
(x0^2 + 1)/x1
1680
(x1^6 + x0^4 + 2*x0^2*x1^2 + 3*x1^4 + 2*x0^2 + 3*x1^2 + 1)/(x0^3*x1^2)
1681
(x1^8 + x0^6 + 2*x0^4*x1^2 + 3*x0^2*x1^4 + 4*x1^6 + 3*x0^4 + 6*x0^2*x1^2 + 6*x1^4 + 3*x0^2 + 4*x1^2 + 1)/(x0^4*x1^3)
1682
(x0^8 + 4*x0^6 + 3*x0^4*x1^2 + 2*x0^2*x1^4 + x1^6 + 6*x0^4 + 6*x0^2*x1^2 + 3*x1^4 + 4*x0^2 + 3*x1^2 + 1)/(x0^3*x1^4)
1683
(x0^6 + 3*x0^4 + 2*x0^2*x1^2 + x1^4 + 3*x0^2 + 2*x1^2 + 1)/(x0^2*x1^3)
1684
"""
1685
mut_iter = self.mutation_class_iter( depth=depth,show_depth=False )
1686
var_class = set()
1687
1688
for seed in mut_iter:
1689
if seed is self:
1690
seed = ClusterSeed(seed)
1691
if not ignore_bipartite_belt and seed.is_bipartite():
1692
bipartition = seed.is_bipartite(return_bipartition=True)
1693
bipartition = (list(bipartition[0]),list(bipartition[1]))
1694
if depth is not infinity:
1695
print "Found a bipartite seed - restarting the depth counter at zero and constructing the variable class using its bipartite belt."
1696
depth_counter = 0
1697
end = False
1698
seed2 = ClusterSeed(seed)
1699
for c in seed._cluster:
1700
if c not in var_class:
1701
yield ClusterVariable( c.parent(), c.numerator(), c.denominator(), mutation_type=self._mutation_type, variable_type='cluster variable' )
1702
var_class = var_class.union( seed._cluster )
1703
1704
init_cluster = set(seed._cluster)
1705
while not end and depth_counter < depth:
1706
depth_counter += 1
1707
seed.mutate(bipartition[0])
1708
seed.mutate(bipartition[1])
1709
if set(seed._cluster) in [set(seed2._cluster),init_cluster]:
1710
end = True
1711
if not end:
1712
for c in seed._cluster:
1713
if c not in var_class:
1714
yield ClusterVariable( c.parent(), c.numerator(), c.denominator(), mutation_type=self._mutation_type, variable_type='cluster variable' )
1715
var_class = var_class.union( seed._cluster )
1716
seed2.mutate(bipartition[1])
1717
seed2.mutate(bipartition[0])
1718
if set(seed2._cluster) in [set(seed._cluster),init_cluster]:
1719
end = True
1720
if not end:
1721
for c in seed2._cluster:
1722
if c not in var_class:
1723
yield ClusterVariable( c.parent(), c.numerator(), c.denominator(), mutation_type=self._mutation_type, variable_type='cluster variable' )
1724
var_class = var_class.union(seed2._cluster)
1725
return
1726
else:
1727
for c in seed._cluster:
1728
if c not in var_class:
1729
yield ClusterVariable( c.parent(), c.numerator(), c.denominator(), mutation_type=self._mutation_type, variable_type='cluster variable' )
1730
var_class = var_class.union(seed._cluster)
1731
1732
def variable_class(self, depth=infinity, ignore_bipartite_belt=False):
1733
r"""
1734
Returns all cluster variables in the mutation class of ``self``.
1735
1736
INPUT:
1737
1738
- ``depth`` -- (default:infinity) integer, only seeds with distance at most depth from self are returned
1739
- ``ignore_bipartite_belt`` -- (default:False) if True, the algorithms does not use the bipartite belt
1740
1741
1742
EXAMPLES:
1743
1744
- for examples see :meth:`variable_class_iter`
1745
1746
TESTS::
1747
1748
sage: A = ClusterSeed(['A',3]).variable_class()
1749
"""
1750
if depth is infinity and not self.is_finite():
1751
raise ValueError('The variable class can - for infinite types - only be computed up to a given depth')
1752
1753
var_iter = self.variable_class_iter( depth=depth, ignore_bipartite_belt=ignore_bipartite_belt )
1754
Vs = [ var for var in var_iter ]
1755
Vs.sort(cmp=cmp)
1756
return Vs
1757
1758
def is_finite( self ):
1759
r"""
1760
Returns True if ``self`` is of finite type.
1761
1762
EXAMPLES::
1763
1764
sage: S = ClusterSeed(['A',3])
1765
sage: S.is_finite()
1766
True
1767
1768
sage: S = ClusterSeed(['A',[2,2],1])
1769
sage: S.is_finite()
1770
False
1771
"""
1772
mt = self.mutation_type()
1773
if type(mt) is str:
1774
return False
1775
else:
1776
return mt.is_finite()
1777
1778
def is_mutation_finite( self, nr_of_checks=None, return_path=False ):
1779
r"""
1780
Returns True if ``self`` is of finite mutation type.
1781
1782
INPUT:
1783
1784
- ``nr_of_checks`` -- (default: None) number of mutations applied. Standard is 500*(number of vertices of self).
1785
- ``return_path`` -- (default: False) if True, in case of self not being mutation finite, a path from self to a quiver with an edge label (a,-b) and a*b > 4 is returned.
1786
1787
ALGORITHM:
1788
1789
- A cluster seed is mutation infinite if and only if every `b_{ij}*b_{ji} > -4`. Thus, we apply random mutations in random directions
1790
1791
WARNING:
1792
1793
- Uses a non-deterministic method by random mutations in various directions.
1794
- In theory, it can return a wrong True.
1795
1796
EXAMPLES::
1797
1798
sage: S = ClusterSeed(['A',10])
1799
sage: S._mutation_type = None
1800
sage: S.is_mutation_finite()
1801
True
1802
1803
sage: S = ClusterSeed([(0,1),(1,2),(2,3),(3,4),(4,5),(5,6),(6,7),(7,8),(2,9)])
1804
sage: S.is_mutation_finite()
1805
False
1806
"""
1807
is_finite, path = is_mutation_finite(copy(self._M),nr_of_checks=nr_of_checks)
1808
if return_path:
1809
return is_finite, path
1810
else:
1811
return is_finite
1812
1813
def mutation_type(self):
1814
r"""
1815
Returns the mutation_type of each connected component of ``self``, if it can be determined.
1816
Otherwise, the mutation type of this component is set to be unknown.
1817
1818
The mutation types of the components are ordered by vertex labels.
1819
1820
WARNING:
1821
1822
- All finite types can be detected,
1823
- All affine types can be detected, EXCEPT affine type D (the algorithm is not yet implemented)
1824
- All exceptional types can be detected.
1825
1826
- Might fail to work if it is used within different Sage processes simultaneously (that happend in the doctesting).
1827
1828
EXAMPLES:
1829
1830
- finite types::
1831
1832
sage: S = ClusterSeed(['A',5])
1833
sage: S._mutation_type = S._quiver._mutation_type = None
1834
sage: S.mutation_type()
1835
['A', 5]
1836
1837
sage: S = ClusterSeed([(0,1),(1,2),(2,3),(3,4)])
1838
sage: S.mutation_type()
1839
['A', 5]
1840
1841
- affine types::
1842
1843
sage: S = ClusterSeed(['E',8,[1,1]]); S
1844
A seed for a cluster algebra of rank 10 of type ['E', 8, [1, 1]]
1845
sage: S._mutation_type = S._quiver._mutation_type = None; S
1846
A seed for a cluster algebra of rank 10
1847
sage: S.mutation_type() # long time
1848
['E', 8, [1, 1]]
1849
1850
- the not yet working affine type D::
1851
1852
sage: S = ClusterSeed(['D',4,1])
1853
sage: S._mutation_type = S._quiver._mutation_type = None
1854
sage: S.mutation_type() # todo: not implemented
1855
['D', 4, 1]
1856
1857
- the exceptional types::
1858
1859
sage: S = ClusterSeed(['X',6])
1860
sage: S._mutation_type = S._quiver._mutation_type = None
1861
sage: S.mutation_type() # long time
1862
['X', 6]
1863
1864
- infinite types::
1865
1866
sage: S = ClusterSeed(['GR',[4,9]])
1867
sage: S._mutation_type = S._quiver._mutation_type = None
1868
sage: S.mutation_type()
1869
'undetermined infinite mutation type'
1870
"""
1871
if self._mutation_type is None:
1872
if self._quiver is None:
1873
self.quiver()
1874
self._mutation_type = self._quiver.mutation_type()
1875
return self._mutation_type
1876
1877
def greedy(self, a1, a2, method='by_recursion'):
1878
r"""
1879
Returns the greedy element `x[a_1,a_2]` assuming that self is rank two.
1880
1881
The third input can be 'by_recursion', 'by_combinatorics', or
1882
'just_numbers' to specify if the user wants the element
1883
computed by the recurrence, combinatorial formula, or wants to
1884
set `x_1` and `x_2` to be one.
1885
1886
See [LeeLiZe]_ for more details.
1887
1888
EXAMPLES::
1889
1890
sage: S = ClusterSeed(['R2', [3, 3]])
1891
sage: S.greedy(4, 4)
1892
(x0^12 + x1^12 + 4*x0^9 + 4*x1^9 + 6*x0^6 + 4*x0^3*x1^3 + 6*x1^6 + 4*x0^3 + 4*x1^3 + 1)/(x0^4*x1^4)
1893
sage: S.greedy(4, 4, 'by_combinatorics')
1894
(x0^12 + x1^12 + 4*x0^9 + 4*x1^9 + 6*x0^6 + 4*x0^3*x1^3 + 6*x1^6 + 4*x0^3 + 4*x1^3 + 1)/(x0^4*x1^4)
1895
sage: S.greedy(4, 4, 'just_numbers')
1896
35
1897
sage: S = ClusterSeed(['R2', [2, 2]])
1898
sage: S.greedy(1, 2)
1899
(x0^4 + 2*x0^2 + x1^2 + 1)/(x0*x1^2)
1900
sage: S.greedy(1, 2, 'by_combinatorics')
1901
(x0^4 + 2*x0^2 + x1^2 + 1)/(x0*x1^2)
1902
1903
REFERENCES:
1904
1905
.. [LeeLiZe] Lee-Li-Zelevinsky, Greedy elements in rank 2
1906
cluster algebras, :arxiv:`1208.2391`
1907
"""
1908
if self.b_matrix().dimensions() == (2, 2):
1909
b = abs(self.b_matrix()[0, 1])
1910
c = abs(self.b_matrix()[1, 0])
1911
if method == 'by_recursion':
1912
ans = self.x(0)**(-a1)*self.x(1)**(-a2)
1913
for p in range(max(a2, 0)+1):
1914
for q in range(max(a1, 0)+1):
1915
if p != 0 or q != 0:
1916
ans += self._R(coeff_recurs(p, q, a1, a2, b, c))*self.x(0)**(b*p-a1)*self.x(1)**(c*q-a2)
1917
return(ans)
1918
elif method == 'by_combinatorics':
1919
if b == 0:
1920
S = ClusterSeed([['A', 1], ['A', 1]])
1921
else:
1922
S = ClusterSeed(['R2', [b, b]])
1923
ans = 0
1924
if a1 >= a2:
1925
PS = PathSubset(a1, a2)
1926
elif a1 < a2:
1927
PS = PathSubset(a2, a1)
1928
from sage.combinat.subset import Subsets
1929
for T in Subsets(PS):
1930
if a1 >= a2:
1931
if is_LeeLiZel_allowable(T, a1, a2, b, c):
1932
oddT = set(T).intersection(PathSubset(a1, 0))
1933
evenT = set(T).symmetric_difference(oddT)
1934
ans = ans + S.x(0)**(b*len(evenT)) * S.x(1)**(c*len(oddT))
1935
elif a1 < a2:
1936
if is_LeeLiZel_allowable(T, a2, a1, b, c):
1937
oddT = set(T).intersection(PathSubset(a2, 0))
1938
evenT = set(T).symmetric_difference(oddT)
1939
ans = ans + S.x(0)**(b*len(oddT)) * S.x(1)**(c*len(evenT))
1940
ans = ans*S.x(0)**(-a1)*S.x(1)**(-a2)
1941
return ans
1942
elif method == 'just_numbers':
1943
ans = 1
1944
for p in range(max(a2, 0)+1):
1945
for q in range(max(a1, 0)+1):
1946
if p != 0 or q != 0:
1947
ans += coeff_recurs(p, q, a1, a2, b, c)
1948
return(ans)
1949
else:
1950
raise ValueError("The third input should be 'by_recursion', "
1951
"'by_combinatorics', or 'just_numbers'.")
1952
else:
1953
raise ValueError("Greedy elements are only currently "
1954
"defined for cluster seeds of rank two.")
1955
1956
def _bino(n, k):
1957
"""
1958
Binomial coefficient which we define as zero for negative n.
1959
1960
EXAMPLES::
1961
1962
sage: from sage.combinat.cluster_algebra_quiver.cluster_seed import _bino
1963
sage: _bino(3, 2)
1964
3
1965
sage: _bino(-3, 2)
1966
0
1967
"""
1968
if n >= 0:
1969
from sage.rings.arith import binomial
1970
return binomial(n, k)
1971
else:
1972
return 0
1973
1974
def coeff_recurs(p, q, a1, a2, b, c):
1975
"""
1976
Coefficients in Laurent expansion of greedy element, as defined by recursion.
1977
1978
EXAMPLES::
1979
1980
sage: from sage.combinat.cluster_algebra_quiver.cluster_seed import coeff_recurs
1981
sage: coeff_recurs(1, 1, 5, 5, 3, 3)
1982
10
1983
"""
1984
if p == 0 and q == 0:
1985
return 1
1986
elif p < 0 or q < 0:
1987
return 0
1988
else:
1989
if c*a1*q <= b*a2*p:
1990
return sum((-1)**(k-1)*coeff_recurs(p-k, q, a1, a2, b, c)*_bino(a2-c*q+k-1, k)
1991
for k in range(1, p+1))
1992
else:
1993
return sum((-1)**(k-1)*coeff_recurs(p, q-k, a1, a2, b, c)*_bino(a1-b*p+k-1, k)
1994
for k in range(1, q+1))
1995
1996
def PathSubset(n,m):
1997
r"""
1998
Encodes a *maximal* Dyck path from (0,0) to (n,m) (for n >= m >= 0) as a subset of {0,1,2,..., 2n-1}.
1999
The encoding is given by indexing horizontal edges by odd numbers and vertical edges by evens.
2000
2001
The horizontal between (i,j) and (i+1,j) is indexed by the odd number 2*i+1.
2002
The vertical between (i,j) and (i,j+1) is indexed by the even number 2*j.
2003
2004
EXAMPLES::
2005
2006
sage: from sage.combinat.cluster_algebra_quiver.cluster_seed import PathSubset
2007
sage: PathSubset(4,0)
2008
set([1, 3, 5, 7])
2009
sage: PathSubset(4,1)
2010
set([1, 3, 5, 6, 7])
2011
sage: PathSubset(4,2)
2012
set([1, 2, 3, 5, 6, 7])
2013
sage: PathSubset(4,3)
2014
set([1, 2, 3, 4, 5, 6, 7])
2015
sage: PathSubset(4,4)
2016
set([0, 1, 2, 3, 4, 5, 6, 7])
2017
"""
2018
from sage.misc.misc import union
2019
from sage.functions.other import floor
2020
S = [ ]
2021
for i in range(n):
2022
S = union(S, [2*i+1])
2023
if m > 0:
2024
for j in range(n):
2025
if floor((j+1)*m/n) - floor(j*m/n) == 1:
2026
S = union(S, [2*j])
2027
return set(S)
2028
2029
def SetToPath(T):
2030
r"""
2031
Rearranges the encoding for a *maximal* Dyck path (as a set) so that it is a list in the proper order of the edges.
2032
2033
EXAMPLES::
2034
2035
sage: from sage.combinat.cluster_algebra_quiver.cluster_seed import PathSubset
2036
sage: from sage.combinat.cluster_algebra_quiver.cluster_seed import SetToPath
2037
sage: SetToPath(PathSubset(4,0))
2038
[1, 3, 5, 7]
2039
sage: SetToPath(PathSubset(4,1))
2040
[1, 3, 5, 7, 6]
2041
sage: SetToPath(PathSubset(4,2))
2042
[1, 3, 2, 5, 7, 6]
2043
sage: SetToPath(PathSubset(4,3))
2044
[1, 3, 2, 5, 4, 7, 6]
2045
sage: SetToPath(PathSubset(4,4))
2046
[1, 0, 3, 2, 5, 4, 7, 6]
2047
"""
2048
n = (max(T)+1)/2
2049
ans = [1]
2050
for i in range(n-1):
2051
if 2*i in T:
2052
ans.append(2*i)
2053
ans.append(2*i+3)
2054
if 2*n-2 in T:
2055
ans.append(2*n-2)
2056
return ans
2057
2058
def is_LeeLiZel_allowable(T,n,m,b,c):
2059
"""
2060
Check if the subset T contributes to the computation of the greedy
2061
element x[m,n] in the rank two (b,c)-cluster algebra.
2062
2063
This uses the conditions of Lee-Li-Zelevinsky's paper [LeeLiZe]_.
2064
2065
EXAMPLES::
2066
2067
sage: from sage.combinat.cluster_algebra_quiver.cluster_seed import is_LeeLiZel_allowable
2068
sage: is_LeeLiZel_allowable({1,3,2,5,7,6},4,2,6,6)
2069
False
2070
sage: is_LeeLiZel_allowable({1,2,5},3,3,1,1)
2071
True
2072
"""
2073
horiz = set(T).intersection( PathSubset(n, 0))
2074
vert = set(T).symmetric_difference(horiz)
2075
if len(horiz) == 0 or len(vert) == 0:
2076
return True
2077
else:
2078
Latt = SetToPath(PathSubset(n, m))
2079
for u in horiz:
2080
from sage.combinat.words.word import Word
2081
from sage.modules.free_module_element import vector
2082
WW = Word(Latt)
2083
LattCycled = vector(WW.conjugate(Latt.index(u))).list()
2084
for v in vert:
2085
uv_okay = False
2086
for A in range(LattCycled.index(v)):
2087
EA = []
2088
AF = copy(LattCycled)
2089
for i in range(LattCycled.index(v), len(LattCycled)-1):
2090
AF.pop()
2091
AF.reverse()
2092
for i in range(A+1):
2093
EA.append(LattCycled[i])
2094
AF.pop()
2095
AF.reverse()
2096
nAF1 = 0
2097
for i in range(len(AF)):
2098
if AF[i] % 2 == 1:
2099
nAF1 += 1
2100
nAF2 = 0
2101
for i in range(len(AF)):
2102
if AF[i] % 2 == 0 and AF[i] in vert:
2103
nAF2 += 1
2104
nEA2 = 0
2105
for i in range(len(EA)):
2106
if EA[i] % 2 == 0:
2107
nEA2 += 1
2108
nEA1 = 0
2109
for i in range(len(EA)):
2110
if EA[i] % 2 == 1 and EA[i] in horiz:
2111
nEA1 += 1
2112
if nAF1 == b*nAF2 or nEA2 == c*nEA1:
2113
uv_okay = True
2114
if uv_okay == False:
2115
return False
2116
return True
2117
2118
2119
class ClusterVariable(FractionFieldElement):
2120
r"""
2121
This class is a thin wrapper for cluster variables in cluster seeds.
2122
2123
It provides the extra feature to store if a variable is frozen or not.
2124
2125
- the associated positive root::
2126
2127
sage: S = ClusterSeed(['A',3])
2128
sage: for T in S.variable_class_iter(): print T, T.almost_positive_root()
2129
x0 -alpha[1]
2130
x1 -alpha[2]
2131
x2 -alpha[3]
2132
(x1 + 1)/x0 alpha[1]
2133
(x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2) alpha[1] + alpha[2] + alpha[3]
2134
(x1 + 1)/x2 alpha[3]
2135
(x0*x2 + x1 + 1)/(x0*x1) alpha[1] + alpha[2]
2136
(x0*x2 + 1)/x1 alpha[2]
2137
(x0*x2 + x1 + 1)/(x1*x2) alpha[2] + alpha[3]
2138
"""
2139
def __init__( self, parent, numerator, denominator, coerce=True, reduce=True, mutation_type=None, variable_type=None ):
2140
r"""
2141
Initializes a cluster variable in the same way that elements in the field of rational functions are initialized.
2142
2143
.. see also:: :class:`Fraction Field of Multivariate Polynomial Ring`
2144
2145
TESTS::
2146
2147
sage: S = ClusterSeed(['A',2])
2148
sage: for f in S.cluster():
2149
... print type(f)
2150
<class 'sage.combinat.cluster_algebra_quiver.cluster_seed.ClusterVariable'>
2151
<class 'sage.combinat.cluster_algebra_quiver.cluster_seed.ClusterVariable'>
2152
2153
sage: S.variable_class()
2154
[(x0 + x1 + 1)/(x0*x1), (x1 + 1)/x0, (x0 + 1)/x1, x1, x0]
2155
"""
2156
FractionFieldElement.__init__( self, parent, numerator, denominator, coerce=coerce, reduce=reduce )
2157
self._mutation_type = mutation_type
2158
self._variable_type = variable_type
2159
2160
def almost_positive_root( self ):
2161
r"""
2162
Returns the *almost positive root* associated to ``self`` if ``self`` is of finite type.
2163
2164
EXAMPLES::
2165
2166
sage: S = ClusterSeed(['A',3])
2167
sage: for T in S.variable_class_iter(): print T, T.almost_positive_root()
2168
x0 -alpha[1]
2169
x1 -alpha[2]
2170
x2 -alpha[3]
2171
(x1 + 1)/x0 alpha[1]
2172
(x1^2 + x0*x2 + 2*x1 + 1)/(x0*x1*x2) alpha[1] + alpha[2] + alpha[3]
2173
(x1 + 1)/x2 alpha[3]
2174
(x0*x2 + x1 + 1)/(x0*x1) alpha[1] + alpha[2]
2175
(x0*x2 + 1)/x1 alpha[2]
2176
(x0*x2 + x1 + 1)/(x1*x2) alpha[2] + alpha[3]
2177
"""
2178
if self._variable_type == 'frozen variable':
2179
raise ValueError('The variable is frozen.')
2180
if type(self._mutation_type) is str:
2181
raise ValueError('The cluster algebra for %s is not of finite type.'%self._repr_())
2182
else:
2183
if self._mutation_type is None:
2184
self._mutation_type = self.parent().mutation_type()
2185
if self._mutation_type.is_finite():
2186
from sage.combinat.root_system.root_system import RootSystem
2187
# the import above is used in the line below
2188
exec "Phi = RootSystem("+self._mutation_type._repr_()+")"
2189
Phiplus = Phi.root_lattice().simple_roots()
2190
if self.denominator() == 1:
2191
return -Phiplus[ self.numerator().degrees().index(1) + 1 ]
2192
else:
2193
root = self.denominator().degrees()
2194
return sum( [ root[i]*Phiplus[ i+1 ] for i in range(len(root)) ] )
2195
else:
2196
raise ValueError('The cluster algebra for %s is not of finite type.'%self._repr_())
2197
2198